/*
* Delphes: a framework for fast simulation of a generic collider experiment
* Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
*
* This program is free software: you can redistribute it and/or modify
* it under the terms of the GNU General Public License as published by
* the Free Software Foundation, either version 3 of the License, or
* (at your option) any later version.
*
* This program is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU General Public License for more details.
*
* You should have received a copy of the GNU General Public License
* along with this program. If not, see .
*/
#include
#include
#include
#include
#include "TROOT.h"
#include "TSystem.h"
#include "TApplication.h"
#include "TString.h"
#include "TH1.h"
#include "TH2.h"
#include "TMath.h"
#include "TStyle.h"
#include "TGraph.h"
#include "TCanvas.h"
#include "THStack.h"
#include "TLegend.h"
#include "TPaveText.h"
#include "TClonesArray.h"
#include "TLorentzVector.h"
#include "TGraphErrors.h"
#include "TMultiGraph.h"
#include "classes/DelphesClasses.h"
#include "ExRootAnalysis/ExRootTreeReader.h"
#include "ExRootAnalysis/ExRootTreeWriter.h"
#include "ExRootAnalysis/ExRootTreeBranch.h"
#include "ExRootAnalysis/ExRootResult.h"
#include "ExRootAnalysis/ExRootUtilities.h"
using namespace std;
//------------------------------------------------------------------------------
static const int Nbins = 50;
int objStyle = 1;
int trackStyle = 7;
int towerStyle = 3;
Color_t objColor = kBlack;
Color_t trackColor = kBlack;
Color_t towerColor = kBlack;
double effLegXmin = 0.22;
double effLegXmax = 0.7;
double effLegYmin = 0.22;
double effLegYmax = 0.5;
double resLegXmin = 0.62;
double resLegXmax = 0.9;
double resLegYmin = 0.52;
double resLegYmax = 0.85;
double topLeftLegXmin = 0.22;
double topLeftLegXmax = 0.7;
double topLeftLegYmin = 0.52;
double topLeftLegYmax = 0.85;
unsigned int k;
struct resolPlot
{
TH1* resolHist;
double ptmin;
double ptmax;
double etamin;
double etamax;
double xmin;
double xmax;
TString obj;
resolPlot();
resolPlot(double ptdown, double ptup, TString object);
resolPlot(double etadown, double etaup, double ptdown, double ptup, TString object);
void set(double ptdown, double ptup, TString object, double xmin = 0, double xmax = 2);
void set(double etadown, double etaup, double ptdown, double ptup, TString object, double xmin = 0, double xmax = 2);
void print(){std::cout << ptmin << std::endl;}
};
resolPlot::resolPlot()
{
}
resolPlot::resolPlot(double ptdown, double ptup, TString object)
{
this->set(ptdown,ptup,object);
}
resolPlot::resolPlot(double etadown, double etaup, double ptdown, double ptup, TString object)
{
this->set(etadown, etaup, ptdown, ptup, object);
}
void resolPlot::set(double ptdown, double ptup, TString object, double xmin, double xmax)
{
ptmin = ptdown;
ptmax = ptup;
obj = object;
resolHist = new TH1D(obj+"_delta_pt_"+Form("%4.2f",ptmin)+"_"+Form("%4.2f",ptmax), obj+"_delta_pt_"+Form("%4.2f",ptmin)+"_"+Form("%4.2f",ptmax), 1000, xmin, xmax);
}
void resolPlot::set(double etadown, double etaup, double ptdown, double ptup, TString object, double xmin, double xmax)
{
etamin = etadown;
etamax = etaup;
ptmin = ptdown;
ptmax = ptup;
obj = object;
resolHist = new TH1D(obj+"_delta_pt_"+Form("%4.2f",ptmin)+"_"+Form("%4.2f",ptmax)+"_"+Form("%4.2f",etamin)+"_"+Form("%4.2f",etamax), obj+"_delta_pt_"+Form("%4.2f",ptmin)+"_"+Form("%4.2f",ptmax)+"_"+Form("%4.2f",etamin)+"_"+Form("%4.2f",etamax), 1000, xmin, xmax);
}
void HistogramsCollection(std::vector *histos, double ptmin, double ptmax, TString obj, double xmin = 0, double xmax = 2)
{
double width;
double ptdown;
double ptup;
resolPlot ptemp;
for (int i = 0; i < Nbins; i++)
{
width = (ptmax - ptmin) / Nbins;
ptdown = TMath::Power(10,ptmin + i * width );
ptup = TMath::Power(10,ptmin + (i+1) * width );
ptemp.set(ptdown, ptup, obj, xmin, xmax);
histos->push_back(ptemp);
}
}
void HistogramsCollectionVsEta(std::vector *histos, double etamin, double etamax, double ptmin, double ptmax, TString obj, double xmin = 0, double xmax = 2)
{
resolPlot ptemp;
double width;
double etadown;
double etaup;
for (int i = 0; i < Nbins; i++)
{
width = (etamax - etamin) / Nbins;
etadown = etamin + i * width;
etaup = etamin + (i+1) * width;
ptemp.set(etadown, etaup, ptmin, ptmax, obj, xmin, xmax);
histos->push_back(ptemp);
}
}
//------------------------------------------------------------------------------
class ExRootResult;
class ExRootTreeReader;
//------------------------------------------------------------------------------
void BinLogX(TH1*h)
{
TAxis *axis = h->GetXaxis();
int bins = axis->GetNbins();
Axis_t from = axis->GetXmin();
Axis_t to = axis->GetXmax();
Axis_t width = (to - from) / bins;
Axis_t *new_bins = new Axis_t[bins + 1];
for (int i = 0; i <= bins; i++)
{
new_bins[i] = TMath::Power(10, from + i * width);
}
axis->Set(bins, new_bins);
delete new_bins;
}
//------------------------------------------------------------------------------
template
TH1D* GetEffPt(TClonesArray *branchReco, TClonesArray *branchParticle, TString name, int pdgID, double ptmin, double ptmax, double etamin, double etamax, ExRootTreeReader *treeReader)
{
cout << "** Computing Efficiency of reconstructing "<< branchReco->GetName() << " induced by " << branchParticle->GetName() << " with PID " << pdgID << endl;
Long64_t allEntries = treeReader->GetEntries();
GenParticle *particle;
T *recoObj;
TLorentzVector recoMomentum, genMomentum, bestRecoMomentum;
Float_t deltaR;
Float_t pt, eta;
Long64_t entry;
Int_t i, j;
TH1D *histGenPt = new TH1D(name+" gen spectra Pt",name+" gen spectra cen", Nbins, TMath::Log10(ptmin), TMath::Log10(ptmax));
TH1D *histRecoPt = new TH1D(name+" reco spectra Pt",name+" reco spectra cen", Nbins, TMath::Log10(ptmin), TMath::Log10(ptmax));
histGenPt->SetDirectory(0);
histRecoPt->SetDirectory(0);
BinLogX(histGenPt);
BinLogX(histRecoPt);
// Loop over all events
for(entry = 0; entry < allEntries; ++entry)
{
// Load selected branches with data from specified event
treeReader->ReadEntry(entry);
// Loop over all generated particle in event
for(i = 0; i < branchParticle->GetEntriesFast(); ++i)
{
particle = (GenParticle*) branchParticle->At(i);
genMomentum = particle->P4();
deltaR = 999;
pt = genMomentum.Pt();
eta = TMath::Abs(genMomentum.Eta());
if(eta > etamax || eta < etamin ) continue;
if (particle->PID == pdgID && genMomentum.Pt() > ptmin && genMomentum.Pt() < ptmax )
//if (TMath::Abs(particle->PID) == pdgID && (particle->Status>20 && particle->Status <30) && genMomentum.Pt() > ptmin && genMomentum.Pt() < ptmax )
{
// Loop over all reco object in event
for(j = 0; j < branchReco->GetEntriesFast(); ++j)
{
recoObj = (T*)branchReco->At(j);
recoMomentum = recoObj->P4();
//if(Momentum.Px() == 0 && genMomentum.Py() == 0) continue;
// take the closest parton candidate
if(TMath::Abs(pdgID) == 5)
{
Jet *jet = (Jet *)recoObj;
if( !(jet->BTag & (1 << 0)) ) continue;
//if(jet->BTag != ) continue;
}
if(TMath::Abs(pdgID) == 4)
{
Jet *jet = (Jet *)recoObj;
if( !(jet->BTag & (1 << 0)) ) continue;
}
if(TMath::Abs(pdgID) == 1)
{
Jet *jet = (Jet *)recoObj;
if( !(jet->BTag & (1 << 0)) ) continue;
}
if(TMath::Abs(pdgID) == 15)
{
Jet *jet = (Jet *)recoObj;
if(jet->TauTag != 1) continue;
}
if(genMomentum.DeltaR(recoMomentum) < deltaR)
{
deltaR = genMomentum.DeltaR(recoMomentum);
bestRecoMomentum = recoMomentum;
}
}
histGenPt->Fill(pt);
if(deltaR < 0.3 && bestRecoMomentum.Pt()> 0.20*pt ) { histRecoPt->Fill(pt);}
}
}
}
histRecoPt->Sumw2();
histGenPt->Sumw2();
histRecoPt->Divide(histGenPt);
histRecoPt->Scale(100.);
return histRecoPt;
}
template
TH1D* GetEffEta(TClonesArray *branchReco, TClonesArray *branchParticle, TString name, int pdgID, double ptmin, double ptmax, double etamin, double etamax, ExRootTreeReader *treeReader)
{
cout << "** Computing Efficiency of reconstructing "<< branchReco->GetName() << " induced by " << branchParticle->GetName() << " with PID " << pdgID << endl;
Long64_t allEntries = treeReader->GetEntries();
GenParticle *particle;
T *recoObj;
TLorentzVector recoMomentum, genMomentum, bestRecoMomentum;
Float_t deltaR;
Float_t pt, eta;
Long64_t entry;
Int_t i, j;
TH1D *histGenEta = new TH1D(name+" gen spectra Eta",name+" gen spectra", Nbins, etamin, etamax);
TH1D *histRecoEta = new TH1D(name+" reco spectra Eta",name+" reco spectra", Nbins, etamin, etamax);
histGenEta->SetDirectory(0);
histRecoEta->SetDirectory(0);
// Loop over all events
for(entry = 0; entry < allEntries; ++entry)
{
// Load selected branches with data from specified event
treeReader->ReadEntry(entry);
// Loop over all generated particle in event
for(i = 0; i < branchParticle->GetEntriesFast(); ++i)
{
particle = (GenParticle*) branchParticle->At(i);
genMomentum = particle->P4();
deltaR = 999;
pt = genMomentum.Pt();
eta = genMomentum.Eta();
if(pt > ptmax || pt < ptmin ) continue;
if (particle->PID == pdgID && genMomentum.Pt() > ptmin && genMomentum.Pt() < ptmax )
//if (TMath::Abs(particle->PID) == pdgID && (particle->Status>20 && particle->Status <30) && genMomentum.Pt() > ptmin && genMomentum.Pt() < ptmax )
{
// Loop over all reco object in event
for(j = 0; j < branchReco->GetEntriesFast(); ++j)
{
recoObj = (T*)branchReco->At(j);
recoMomentum = recoObj->P4();
// this is simply to avoid warnings from initial state particle
// having infite rapidity ...
//if(Momentum.Px() == 0 && genMomentum.Py() == 0) continue;
// take the closest parton candidate
if(TMath::Abs(pdgID) == 5)
{
Jet *jet = (Jet *)recoObj;
if( !(jet->BTag & (1 << 0)) ) continue;
}
if(TMath::Abs(pdgID) == 4)
{
Jet *jet = (Jet *)recoObj;
if( !(jet->BTag & (1 << 0)) ) continue;
}
if(TMath::Abs(pdgID) == 1)
{
Jet *jet = (Jet *)recoObj;
if( !(jet->BTag & (1 << 0)) ) continue;
}
if(TMath::Abs(pdgID) == 15)
{
Jet *jet = (Jet *)recoObj;
if(jet->TauTag != 1) continue;
}
if(genMomentum.DeltaR(recoMomentum) < deltaR)
{
deltaR = genMomentum.DeltaR(recoMomentum);
bestRecoMomentum = recoMomentum;
}
}
histGenEta->Fill(eta);
if(deltaR < 0.3) { histRecoEta->Fill(eta); }
}
}
}
histRecoEta->Sumw2();
histGenEta->Sumw2();
histRecoEta->Divide(histGenEta);
histRecoEta->Scale(100.);
return histRecoEta;
}
//------------------------------------------------------------------------------
template
TH1D* GetJetEffPt(TClonesArray *branchJet, TString name, int pdgID, double ptmin, double ptmax, double etamin, double etamax, ExRootTreeReader *treeReader)
{
cout << "** Computing Efficiency of reconstructing "<< branchJet->GetName() << " with PID " << pdgID << endl;
Long64_t allEntries = treeReader->GetEntries();
T *recoObj;
TLorentzVector recoMomentum, genMomentum, bestRecoMomentum;
Float_t pt, eta;
Long64_t entry;
Int_t j;
TH1D *histGenPt = new TH1D(name+" gen spectra Pt",name+" gen spectra cen", Nbins, TMath::Log10(ptmin), TMath::Log10(ptmax));
TH1D *histRecoPt = new TH1D(name+" reco spectra Pt",name+" reco spectra cen", Nbins, TMath::Log10(ptmin), TMath::Log10(ptmax));
histGenPt->SetDirectory(0);
histRecoPt->SetDirectory(0);
BinLogX(histGenPt);
BinLogX(histRecoPt);
// Loop over all events
for(entry = 0; entry < allEntries; ++entry)
{
// Load selected branches with data from specified event
treeReader->ReadEntry(entry);
// Loop over all reco object in event
for(j = 0; j < branchJet->GetEntriesFast(); ++j)
{
recoObj = (T*)branchJet->At(j);
recoMomentum = recoObj->P4();
pt = recoMomentum.Pt();
eta = TMath::Abs(recoMomentum.Eta());
Jet *jet = (Jet *)recoObj;
if(eta > etamax || eta < etamin ) continue;
if(pt < ptmin || pt > ptmax) continue;
Int_t flavor = jet->Flavor;
if(flavor == 21) flavor = 0;
if(TMath::Abs(pdgID) == 1)
{
if(flavor < 4)
{
histGenPt->Fill(pt);
if( jet->BTag & (1 << 0) ) histRecoPt->Fill(pt);
}
}
if(TMath::Abs(pdgID) == 4)
{
if(flavor == 4)
{
histGenPt->Fill(pt);
if( jet->BTag & (1 << 0) ) histRecoPt->Fill(pt);
}
}
if(TMath::Abs(pdgID) == 5)
{
if(flavor == 5)
{
histGenPt->Fill(pt);
if( jet->BTag & (1 << 0) ) histRecoPt->Fill(pt);
}
}
}
}
histRecoPt->Sumw2();
histGenPt->Sumw2();
histRecoPt->Divide(histGenPt);
histRecoPt->Scale(100.);
return histRecoPt;
}
// ------------------------------------------------------------------------------------------------------------------------------------------------------
template
TH1D* GetJetEffEta(TClonesArray *branchJet, TString name, int pdgID, double ptmin, double ptmax, double etamin, double etamax, ExRootTreeReader *treeReader)
{
cout << "** Computing Efficiency of reconstructing "<< branchJet->GetName() << " with PID " << pdgID << endl;
Long64_t allEntries = treeReader->GetEntries();
T *recoObj;
TLorentzVector recoMomentum, genMomentum, bestRecoMomentum;
Float_t pt, eta;
Long64_t entry;
Int_t j;
TH1D *histGenEta = new TH1D(name+" gen spectra Eta",name+" gen spectra", Nbins, etamin, etamax);
TH1D *histRecoEta = new TH1D(name+" reco spectra Eta",name+" reco spectra", Nbins, etamin, etamax);
histGenEta->SetDirectory(0);
histRecoEta->SetDirectory(0);
// Loop over all events
for(entry = 0; entry < allEntries; ++entry)
{
// Load selected branches with data from specified event
treeReader->ReadEntry(entry);
// Loop over all reco object in event
for(j = 0; j < branchJet->GetEntriesFast(); ++j)
{
recoObj = (T*)branchJet->At(j);
recoMomentum = recoObj->P4();
pt = recoMomentum.Pt();
eta = recoMomentum.Eta();
Jet *jet = (Jet *)recoObj;
if(eta > etamax || eta < etamin ) continue;
if(pt < ptmin || pt > ptmax) continue;
Int_t flavor = jet->Flavor;
if(flavor == 21) flavor = 0;
if(TMath::Abs(pdgID) == 1)
{
if(flavor == 1 || flavor == 21)
{
histGenEta->Fill(eta);
if( jet->BTag & (1 << 0) ) histRecoEta->Fill(eta);
}
}
if(TMath::Abs(pdgID) == 4)
{
if(flavor == 4)
{
histGenEta->Fill(eta);
if( jet->BTag & (1 << 0) ) histRecoEta->Fill(eta);
}
}
if(TMath::Abs(pdgID) == 5)
{
if(flavor == 5)
{
histGenEta->Fill(eta);
if( jet->BTag & (1 << 0) ) histRecoEta->Fill(eta);
}
}
}
}
histRecoEta->Sumw2();
histGenEta->Sumw2();
histRecoEta->Divide(histGenEta);
histRecoEta->Scale(100.);
return histRecoEta;
}
// -----------------------------------------------------------------------------------------------------------------------------------------------------
template
TH1D* GetTauEffPt(TClonesArray *branchReco, TClonesArray *branchParticle, TString name, int pdgID, double ptmin, double ptmax, double etamin, double etamax, ExRootTreeReader *treeReader)
{
cout << "** Computing Efficiency of reconstructing "<< branchReco->GetName() << " induced by " << branchParticle->GetName() << " with PID " << pdgID << endl;
Long64_t allEntries = treeReader->GetEntries();
GenParticle *particle;
T *recoObj;
TLorentzVector recoMomentum, genMomentum, bestRecoMomentum;
Float_t deltaR;
Float_t pt, eta;
Long64_t entry;
Int_t i, j;
TH1D *histGenPt = new TH1D(name+" gen spectra Pt",name+" gen spectra cen", Nbins, TMath::Log10(ptmin), TMath::Log10(ptmax));
TH1D *histRecoPt = new TH1D(name+" reco spectra Pt",name+" reco spectra cen", Nbins, TMath::Log10(ptmin), TMath::Log10(ptmax));
histGenPt->SetDirectory(0);
histRecoPt->SetDirectory(0);
BinLogX(histGenPt);
BinLogX(histRecoPt);
// Loop over all events
for(entry = 0; entry < allEntries; ++entry)
{
// Load selected branches with data from specified event
treeReader->ReadEntry(entry);
// Loop over all generated particle in event
for(i = 0; i < branchParticle->GetEntriesFast(); ++i)
{
particle = (GenParticle*) branchParticle->At(i);
genMomentum = particle->P4();
deltaR = 999;
pt = genMomentum.Pt();
eta = TMath::Abs(genMomentum.Eta());
if(eta > etamax || eta < etamin ) continue;
if (particle->PID == pdgID && genMomentum.Pt() > ptmin && genMomentum.Pt() < ptmax )
{
// Loop over all reco object in event
for(j = 0; j < branchReco->GetEntriesFast(); ++j)
{
recoObj = (T*)branchReco->At(j);
recoMomentum = recoObj->P4();
// this is simply to avoid warnings from initial state particle
// having infite rapidity ...
//if(Momentum.Px() == 0 && genMomentum.Py() == 0) continue;
if(TMath::Abs(pdgID) == 1)
{
Jet *jet = (Jet *)recoObj;
if( jet->TauTag != 1 ) continue;
}
if(TMath::Abs(pdgID) == 15)
{
Jet *jet = (Jet *)recoObj;
if(jet->TauTag != 1) continue;
}
if(genMomentum.DeltaR(recoMomentum) < deltaR)
{
deltaR = genMomentum.DeltaR(recoMomentum);
bestRecoMomentum = recoMomentum;
}
}
histGenPt->Fill(pt);
if(deltaR < 0.3) { histRecoPt->Fill(pt); }
}
}
}
histRecoPt->Sumw2();
histGenPt->Sumw2();
histRecoPt->Divide(histGenPt);
histRecoPt->Scale(100.);
if(TMath::Abs(pdgID) == 15) histRecoPt->Scale(1/0.648);
return histRecoPt;
}
template
TH1D* GetTauEffEta(TClonesArray *branchReco, TClonesArray *branchParticle, TString name, int pdgID, double ptmin, double ptmax, double etamin, double etamax, ExRootTreeReader *treeReader)
{
cout << "** Computing Efficiency of reconstructing "<< branchReco->GetName() << " induced by " << branchParticle->GetName() << " with PID " << pdgID << endl;
Long64_t allEntries = treeReader->GetEntries();
GenParticle *particle;
T *recoObj;
TLorentzVector recoMomentum, genMomentum, bestRecoMomentum;
Float_t deltaR;
Float_t pt, eta;
Long64_t entry;
Int_t i, j;
TH1D *histGenEta = new TH1D(name+" gen spectra Eta",name+" gen spectra", Nbins, etamin, etamax);
TH1D *histRecoEta = new TH1D(name+" reco spectra Eta",name+" reco spectra", Nbins, etamin, etamax);
histGenEta->SetDirectory(0);
histRecoEta->SetDirectory(0);
// Loop over all events
for(entry = 0; entry < allEntries; ++entry)
{
// Load selected branches with data from specified event
treeReader->ReadEntry(entry);
// Loop over all generated particle in event
for(i = 0; i < branchParticle->GetEntriesFast(); ++i)
{
particle = (GenParticle*) branchParticle->At(i);
genMomentum = particle->P4();
deltaR = 999;
pt = genMomentum.Pt();
eta = genMomentum.Eta();
if(pt > ptmax || pt < ptmin ) continue;
if (particle->PID == pdgID && genMomentum.Pt() > ptmin && genMomentum.Pt() < ptmax )
{
// Loop over all reco object in event
for(j = 0; j < branchReco->GetEntriesFast(); ++j)
{
recoObj = (T*)branchReco->At(j);
recoMomentum = recoObj->P4();
// this is simply to avoid warnings from initial state particle
// having infite rapidity ...
//if(Momentum.Px() == 0 && genMomentum.Py() == 0) continue;
if(TMath::Abs(pdgID) == 1)
{
Jet *jet = (Jet *)recoObj;
if( jet->TauTag != 1 ) continue;
}
if(TMath::Abs(pdgID) == 15)
{
Jet *jet = (Jet *)recoObj;
if(jet->TauTag != 1) continue;
}
if(genMomentum.DeltaR(recoMomentum) < deltaR)
{
deltaR = genMomentum.DeltaR(recoMomentum);
bestRecoMomentum = recoMomentum;
}
}
histGenEta->Fill(eta);
if(deltaR < 0.3) { histRecoEta->Fill(eta); }
}
}
}
histRecoEta->Sumw2();
histGenEta->Sumw2();
histRecoEta->Divide(histGenEta);
histRecoEta->Scale(100.);
if(TMath::Abs(pdgID) == 15) histRecoEta->Scale(1/0.648);
return histRecoEta;
}
template
void GetPtres(std::vector *histos, TClonesArray *branchReco, TClonesArray *branchParticle, int pdgID, Double_t etaMin, Double_t etaMax, ExRootTreeReader *treeReader)
{
Long64_t allEntries = treeReader->GetEntries();
cout << "** Computing pt resolution of " << branchReco->GetName() << " induced by " << branchParticle->GetName() << " with PID " << pdgID << endl;
GenParticle *particle;
T* recoObj;
TLorentzVector recoMomentum, genMomentum, bestGenMomentum;
Float_t deltaR;
Float_t pt, eta;
Long64_t entry;
Int_t i, j, bin;
// Loop over all events
for(entry = 0; entry < allEntries; ++entry)
{
// Load selected branches with data from specified event
treeReader->ReadEntry(entry);
// Loop over all reconstructed jets in event
for(i = 0; i < branchReco->GetEntriesFast(); ++i)
{
recoObj = (T*) branchReco->At(i);
recoMomentum = recoObj->P4();
deltaR = 999;
// Loop over all hard partons in event
for(j = 0; j < branchParticle->GetEntriesFast(); ++j)
{
particle = (GenParticle*) branchParticle->At(j);
if (particle->PID == pdgID && particle->Status == 1)
{
genMomentum = particle->P4();
// this is simply to avoid warnings from initial state particle
// having infite rapidity ...
if(genMomentum.Px() == 0 && genMomentum.Py() == 0) continue;
// take the closest parton candidate
if(genMomentum.DeltaR(recoMomentum) < deltaR)
{
deltaR = genMomentum.DeltaR(recoMomentum);
bestGenMomentum = genMomentum;
}
}
}
if(deltaR < 0.3)
{
pt = bestGenMomentum.Pt();
eta = TMath::Abs(bestGenMomentum.Eta());
for (bin = 0; bin < Nbins; bin++)
{
if(pt > histos->at(bin).ptmin && pt < histos->at(bin).ptmax && eta > etaMin && eta < etaMax)
{
histos->at(bin).resolHist->Fill(recoMomentum.Pt()/bestGenMomentum.Pt());
}
}
}
}
}
}
template
void GetEres(std::vector *histos, TClonesArray *branchReco, TClonesArray *branchParticle, int pdgID, Double_t etaMin, Double_t etaMax, ExRootTreeReader *treeReader)
{
Long64_t allEntries = treeReader->GetEntries();
cout << "** Computing e resolution of " << branchReco->GetName() << " induced by " << branchParticle->GetName() << " with PID " << pdgID << endl;
GenParticle *particle;
T* recoObj;
TLorentzVector recoMomentum, genMomentum, bestGenMomentum;
Float_t deltaR;
Float_t e, eta;
Long64_t entry;
Int_t i, j, bin;
// Loop over all events
for(entry = 0; entry < allEntries; ++entry)
{
// Load selected branches with data from specified event
treeReader->ReadEntry(entry);
// Loop over all reconstructed jets in event
for(i = 0; i < branchReco->GetEntriesFast(); ++i)
{
recoObj = (T*) branchReco->At(i);
recoMomentum = recoObj->P4();
deltaR = 999;
// Loop over all hard partons in event
for(j = 0; j < branchParticle->GetEntriesFast(); ++j)
{
particle = (GenParticle*) branchParticle->At(j);
if (particle->PID == pdgID && particle->Status == 1)
{
genMomentum = particle->P4();
// this is simply to avoid warnings from initial state particle
// having infite rapidity ...
if(genMomentum.Px() == 0 && genMomentum.Py() == 0) continue;
// take the closest parton candidate
if(genMomentum.DeltaR(recoMomentum) < deltaR)
{
deltaR = genMomentum.DeltaR(recoMomentum);
bestGenMomentum = genMomentum;
}
}
}
if(deltaR < 0.3)
{
e = bestGenMomentum.E();
eta = TMath::Abs(bestGenMomentum.Eta());
for (bin = 0; bin < Nbins; bin++)
{
if(e > histos->at(bin).ptmin && e < histos->at(bin).ptmax && eta > etaMin && eta < etaMax)
{
histos->at(bin).resolHist->Fill(recoMomentum.E()/bestGenMomentum.E());
}
}
}
}
}
}
template
void GetPtresVsEta(std::vector *histos, TClonesArray *branchReco, TClonesArray *branchParticle, int pdgID, Double_t ptMin, Double_t ptMax, ExRootTreeReader *treeReader)
{
Long64_t allEntries = treeReader->GetEntries();
cout << "** Computing pt resolution of " << branchReco->GetName() << " induced by " << branchParticle->GetName() << " with PID " << pdgID << endl;
GenParticle *particle;
T* recoObj;
TLorentzVector recoMomentum, genMomentum, bestGenMomentum;
Float_t deltaR;
Float_t pt, eta;
Long64_t entry;
Int_t i, j, bin;
// Loop over all events
for(entry = 0; entry < allEntries; ++entry)
{
// Load selected branches with data from specified event
treeReader->ReadEntry(entry);
// Loop over all reconstructed jets in event
for(i = 0; i < branchReco->GetEntriesFast(); ++i)
{
recoObj = (T*) branchReco->At(i);
recoMomentum = recoObj->P4();
deltaR = 999;
// Loop over all hard partons in event
for(j = 0; j < branchParticle->GetEntriesFast(); ++j)
{
particle = (GenParticle*) branchParticle->At(j);
if (particle->PID == pdgID && particle->Status == 1)
{
genMomentum = particle->P4();
// this is simply to avoid warnings from initial state particle
// having infite rapidity ...
if(genMomentum.Px() == 0 && genMomentum.Py() == 0) continue;
// take the closest parton candidate
if(genMomentum.DeltaR(recoMomentum) < deltaR)
{
deltaR = genMomentum.DeltaR(recoMomentum);
bestGenMomentum = genMomentum;
}
}
}
if(deltaR < 0.3)
{
pt = bestGenMomentum.Pt();
eta = bestGenMomentum.Eta();
for (bin = 0; bin < Nbins; bin++)
{
if(eta > histos->at(bin).etamin && eta < histos->at(bin).etamax && pt > ptMin && pt < ptMax)
{
histos->at(bin).resolHist->Fill(recoMomentum.Pt()/bestGenMomentum.Pt());
}
}
}
}
}
}
template
void GetEresVsEta(std::vector *histos, TClonesArray *branchReco, TClonesArray *branchParticle, int pdgID, Double_t eMin, Double_t eMax, ExRootTreeReader *treeReader)
{
Long64_t allEntries = treeReader->GetEntries();
cout << "** Computing E resolution of " << branchReco->GetName() << " induced by " << branchParticle->GetName() << " with PID " << pdgID << endl;
GenParticle *particle;
T* recoObj;
TLorentzVector recoMomentum, genMomentum, bestGenMomentum;
Float_t deltaR;
Float_t e, eta;
Long64_t entry;
Int_t i, j, bin;
// Loop over all events
for(entry = 0; entry < allEntries; ++entry)
{
// Load selected branches with data from specified event
treeReader->ReadEntry(entry);
// Loop over all reconstructed jets in event
for(i = 0; i < branchReco->GetEntriesFast(); ++i)
{
recoObj = (T*) branchReco->At(i);
recoMomentum = recoObj->P4();
deltaR = 999;
// Loop over all hard partons in event
for(j = 0; j < branchParticle->GetEntriesFast(); ++j)
{
particle = (GenParticle*) branchParticle->At(j);
if (particle->PID == pdgID && particle->Status == 1)
{
genMomentum = particle->P4();
// this is simply to avoid warnings from initial state particle
// having infite rapidity ...
if(genMomentum.Px() == 0 && genMomentum.Py() == 0) continue;
// take the closest parton candidate
if(genMomentum.DeltaR(recoMomentum) < deltaR)
{
deltaR = genMomentum.DeltaR(recoMomentum);
bestGenMomentum = genMomentum;
}
}
}
if(deltaR < 0.3)
{
e = bestGenMomentum.E();
eta = bestGenMomentum.Eta();
for (bin = 0; bin < Nbins; bin++)
{
if(eta > histos->at(bin).etamin && eta < histos->at(bin).etamax && e > eMin && e < eMax)
{
histos->at(bin).resolHist->Fill(recoMomentum.E()/bestGenMomentum.E());
}
}
}
}
}
}
void GetJetsEres(std::vector *histos, TClonesArray *branchJet, TClonesArray *branchGenJet, ExRootTreeReader *treeReader, Double_t etaMin, Double_t etaMax)
{
Long64_t allEntries = treeReader->GetEntries();
cout << "** Computing resolution of " << branchJet->GetName() << " induced by " << branchGenJet->GetName() << endl;
Jet *jet, *genjet;
TLorentzVector jetMomentum, genJetMomentum, bestGenJetMomentum;
Float_t deltaR;
Float_t pt, eta;
Long64_t entry;
Int_t i, j, bin;
// Loop over all events
for(entry = 0; entry < allEntries; ++entry)
{
// Load selected branches with data from specified event
treeReader->ReadEntry(entry);
if(entry%10000 == 0) cout << "Event number: "<< entry <GetEntriesFast()); ++i) //branchJet->GetEntriesFast(); ++i)
{
jet = (Jet*) branchJet->At(i);
jetMomentum = jet->P4();
deltaR = 999;
// Loop over all hard partons in event
for(j = 0; j < TMath::Min(2,branchGenJet->GetEntriesFast()); ++j)
{
genjet = (Jet*) branchGenJet->At(j);
genJetMomentum = genjet->P4();
// this is simply to avoid warnings from initial state particle
// having infite rapidity ...
if(genJetMomentum.Px() == 0 && genJetMomentum.Py() == 0) continue;
// take the closest parton candidate
if(genJetMomentum.DeltaR(jetMomentum) < deltaR)
{
deltaR = genJetMomentum.DeltaR(jetMomentum);
bestGenJetMomentum = genJetMomentum;
}
}
if(deltaR < 0.3)
{
pt = genJetMomentum.E();
eta = genJetMomentum.Eta();
for (bin = 0; bin < Nbins; bin++)
{
if(pt > histos->at(bin).ptmin && pt < histos->at(bin).ptmax && eta < etaMax && eta > etaMin)
{
histos->at(bin).resolHist->Fill(jetMomentum.E()/bestGenJetMomentum.E());
}
}
}
}
}
}
void GetJetsEresVsEta(std::vector *histos, TClonesArray *branchJet, TClonesArray *branchGenJet, ExRootTreeReader *treeReader, Double_t eMin, Double_t eMax)
{
Long64_t allEntries = treeReader->GetEntries();
cout << "** Computing resolution of " << branchJet->GetName() << " induced by " << branchGenJet->GetName() << endl;
Jet *jet, *genjet;
TLorentzVector jetMomentum, genJetMomentum, bestGenJetMomentum;
Float_t deltaR;
Float_t pt, eta;
Long64_t entry;
Int_t i, j, bin;
// Loop over all events
for(entry = 0; entry < allEntries; ++entry)
{
// Load selected branches with data from specified event
treeReader->ReadEntry(entry);
if(entry%10000 == 0) cout << "Event number: "<< entry <GetEntriesFast()); ++i) //branchJet->GetEntriesFast(); ++i)
{
jet = (Jet*) branchJet->At(i);
jetMomentum = jet->P4();
deltaR = 999;
// Loop over all hard partons in event
for(j = 0; j < TMath::Min(2,branchGenJet->GetEntriesFast()); ++j)
{
genjet = (Jet*) branchGenJet->At(j);
genJetMomentum = genjet->P4();
// this is simply to avoid warnings from initial state particle
// having infite rapidity ...
if(genJetMomentum.Px() == 0 && genJetMomentum.Py() == 0) continue;
// take the closest parton candidate
if(genJetMomentum.DeltaR(jetMomentum) < deltaR)
{
deltaR = genJetMomentum.DeltaR(jetMomentum);
bestGenJetMomentum = genJetMomentum;
}
}
if(deltaR < 0.3)
{
pt = genJetMomentum.E();
eta = genJetMomentum.Eta();
for (bin = 0; bin < Nbins; bin++)
{
if(eta > histos->at(bin).etamin && eta < histos->at(bin).etamax && pt < eMax && pt > eMin)
{
histos->at(bin).resolHist->Fill(jetMomentum.E()/bestGenJetMomentum.E());
}
}
}
}
}
}
std::pair GausFit(TH1* hist)
{
TF1 *f1 = new TF1("f1", "gaus", hist->GetMean()-2*hist->GetRMS(), hist->GetMean()+2*hist->GetRMS());
hist->Fit("f1","RQ");
TF1 *f2 = new TF1("f2", "gaus", f1->GetParameter(1) - 2*f1->GetParameter(2), f1->GetParameter(1) + 2*f1->GetParameter(2));
hist->Fit("f2","RQ");
Double_t sig = f2->GetParameter(2);
Double_t sigErr = f2->GetParError(2);
delete f1;
delete f2;
return make_pair (sig, sigErr);
}
TGraphErrors EresGraph(std::vector *histos, bool rms = false)
{
Int_t bin;
Int_t count = 0;
TGraphErrors gr = TGraphErrors(Nbins/2);
double val, error;
for (bin = 0; bin < Nbins; bin++)
{
std::pair sigvalues = GausFit(histos->at(bin).resolHist);
if (rms == true)
{
gr.SetPoint(count,(histos->at(bin).ptmin+histos->at(bin).ptmax)/2.0, 100*histos->at(bin).resolHist->GetRMS());
//gr.SetPointError(count,0, 100*sigvalues.second); // to correct
error = 100*histos->at(bin).resolHist->GetRMSError();
val = 100*histos->at(bin).resolHist->GetRMS();
if(error > 0.2*val) error = 0.2*val;
gr.SetPointError(count,0, error); // to correct
}
else
{
gr.SetPoint(count,(histos->at(bin).ptmin+histos->at(bin).ptmax)/2.0, 100*sigvalues.first);
error = 100*sigvalues.second;
val = 100*sigvalues.first;
if(error > 0.2*val) error = 0.2*val;
gr.SetPointError(count,0, error); // to correct
//gr.SetPointError(count,0, 100*sigvalues.second);
}
count++;
}
return gr;
}
TGraphErrors MetResGraph(std::vector *histos, bool rms = false)
{
Int_t bin;
Int_t count = 0;
TGraphErrors gr = TGraphErrors(Nbins/2);
double val, error;
for (bin = 0; bin < Nbins; bin++)
{
std::pair sigvalues = GausFit(histos->at(bin).resolHist);
if (rms == true)
{
gr.SetPoint(count,(histos->at(bin).ptmin+histos->at(bin).ptmax)/2.0, histos->at(bin).resolHist->GetRMS());
error = histos->at(bin).resolHist->GetRMSError();
val = histos->at(bin).resolHist->GetRMS();
if(error > 0.2*val) error = 0.2*val;
gr.SetPointError(count,0, error); // to correct
}
else
{
gr.SetPoint(count,(histos->at(bin).ptmin+histos->at(bin).ptmax)/2.0, sigvalues.first);
val = sigvalues.first;
error = sigvalues.second;
if(error > 0.2*val) error = 0.2*val;
gr.SetPointError(count,0, error);
//gr.SetPointError(count,0, 100*sigvalues.second);
}
count++;
}
return gr;
}
TGraphErrors EresGraphVsEta(std::vector *histos, bool rms = false)
{
Int_t bin;
Int_t count = 0;
TGraphErrors gr = TGraphErrors(Nbins/2);
double val, error;
for (bin = 0; bin < Nbins; bin++)
{
std::pair sigvalues = GausFit(histos->at(bin).resolHist);
if (rms == true)
{
gr.SetPoint(count,(histos->at(bin).etamin+histos->at(bin).etamax)/2.0, histos->at(bin).resolHist->GetRMS());
error = 100*histos->at(bin).resolHist->GetRMSError();
val = 100*histos->at(bin).resolHist->GetRMS();
if(error > 0.2*val) error = 0.2*val;
gr.SetPointError(count,0, error); // to correct
}
else
{
gr.SetPoint(count,(histos->at(bin).etamin+histos->at(bin).etamax)/2.0, 100*sigvalues.first);
val = 100*sigvalues.first;
error = 100*sigvalues.second;
if(error > 0.2*val) error = 0.2*val;
gr.SetPointError(count,0, error);
//gr.SetPointError(count,0, 100*sigvalues.second);
}
count++;
}
return gr;
}
void GetMetres(std::vector *histos, TClonesArray *branchScalarHT, TClonesArray *branchMet, TClonesArray *branchJet, ExRootTreeReader *treeReader)
{
Long64_t allEntries = treeReader->GetEntries();
cout << "** Computing resolution of " << branchMet->GetName() << " vs " << branchScalarHT->GetName() << endl;
MissingET *met;
ScalarHT *scalarHT;
Long64_t entry;
Int_t bin;
Double_t ht;
Jet *jet;
TLorentzVector p1, p2;
// Loop over all events
for(entry = 0; entry < allEntries; ++entry)
{
// Load selected branches with data from specified event
treeReader->ReadEntry(entry);
if(entry%10000 == 0) cout << "Event number: "<< entry <GetEntriesFast() > 1)
{
jet = (Jet*) branchJet->At(0);
p1 = jet->P4();
jet = (Jet*) branchJet->At(1);
p2 = jet->P4();
met = (MissingET*) branchMet->At(0);
scalarHT = (ScalarHT*) branchScalarHT->At(0);
ht = scalarHT->HT;
if(p1.Pt() < 0.75*ht/2) continue;
if(p2.Pt() < 0.75*ht/2) continue;
for (bin = 0; bin < Nbins; bin++)
{
if(ht > histos->at(bin).ptmin && ht < histos->at(bin).ptmax )
{
histos->at(bin).resolHist->Fill(met->P4().Px());
}
}
}
}
}
//------------------------------------------------------------------------------
void addResoGraph(TMultiGraph *mg, TGraphErrors *gr, TLegend *leg, int style, Color_t color, TString text)
{
gr->SetLineWidth(2);
gr->SetLineColor(color);
gr->SetMarkerStyle(style);
gr->SetMarkerColor(color);
gr->SetMarkerSize(1.5);
std::cout << "Adding " << gr->GetName() << std::endl;
mg->Add(gr);
leg->AddEntry(gr,text,"p");
}
void DrawAxis(TMultiGraph *mg, TLegend *leg, double xmin, double xmax, double ymin, double ymax, TString tx, TString ty, bool logx = 0, bool logy = 0)
{
mg->SetMinimum(ymin);
mg->SetMaximum(ymax);
mg->GetXaxis()->SetLimits(xmin,xmax);
mg->GetXaxis()->SetTitle(tx);
mg->GetYaxis()->SetTitle(ty);
mg->GetYaxis()->SetTitleSize(0.07);
mg->GetXaxis()->SetTitleSize(0.07);
mg->GetYaxis()->SetLabelSize(0.06);
mg->GetXaxis()->SetLabelSize(0.06);
mg->GetYaxis()->SetLabelOffset(0.03);
mg->GetYaxis()->SetTitleOffset(1.4);
mg->GetXaxis()->SetTitleOffset(1.4);
mg->GetXaxis()->SetTitleFont(132);
mg->GetYaxis()->SetTitleFont(132);
mg->GetXaxis()->SetLabelFont(132);
mg->GetYaxis()->SetLabelFont(132);
mg->GetYaxis()->SetNdivisions(505);
leg->SetTextFont(132);
leg->SetBorderSize(0);
leg->SetShadowColor(0);
leg->SetFillColor(0);
leg->SetFillStyle(0);
gStyle->SetOptTitle(0);
if(logx) gPad->SetLogx();
if(logy) gPad->SetLogy();
//gPad->SetGridx();
//gPad->SetGridy();
gPad->SetBottomMargin(0.2);
gPad->SetLeftMargin(0.2);
gPad->Modified();
gPad->Update();
}
void Validation(const char *inputFilePion,
const char *inputFileElectron,
const char *inputFileMuon,
const char *inputFilePhoton,
const char *inputFileNeutralHadron,
const char *inputFileJet,
const char *inputFileBJet,
const char *inputFileCJet,
const char *inputFileTauJet,
const char *outputFile,
const char *version)
{
TChain *chainPion = new TChain("Delphes");
chainPion->Add(inputFilePion);
ExRootTreeReader *treeReaderPion = new ExRootTreeReader(chainPion);
TChain *chainElectron = new TChain("Delphes");
chainElectron->Add(inputFileElectron);
ExRootTreeReader *treeReaderElectron = new ExRootTreeReader(chainElectron);
TChain *chainMuon = new TChain("Delphes");
chainMuon->Add(inputFileMuon);
ExRootTreeReader *treeReaderMuon = new ExRootTreeReader(chainMuon);
TChain *chainPhoton = new TChain("Delphes");
chainPhoton->Add(inputFilePhoton);
ExRootTreeReader *treeReaderPhoton = new ExRootTreeReader(chainPhoton);
TChain *chainNeutralHadron = new TChain("Delphes");
chainNeutralHadron->Add(inputFileNeutralHadron);
ExRootTreeReader *treeReaderNeutralHadron = new ExRootTreeReader(chainNeutralHadron);
TChain *chainJet = new TChain("Delphes");
chainJet->Add(inputFileJet);
ExRootTreeReader *treeReaderJet = new ExRootTreeReader(chainJet);
TChain *chainBJet = new TChain("Delphes");
chainBJet->Add(inputFileBJet);
ExRootTreeReader *treeReaderBJet = new ExRootTreeReader(chainBJet);
TChain *chainCJet = new TChain("Delphes");
chainCJet->Add(inputFileCJet);
ExRootTreeReader *treeReaderCJet = new ExRootTreeReader(chainCJet);
TChain *chainTauJet = new TChain("Delphes");
chainTauJet->Add(inputFileTauJet);
ExRootTreeReader *treeReaderTauJet = new ExRootTreeReader(chainTauJet);
TClonesArray *branchParticleElectron = treeReaderElectron->UseBranch("Particle");
TClonesArray *branchTrackElectron = treeReaderElectron->UseBranch("Track");
TClonesArray *branchElectron = treeReaderElectron->UseBranch("Electron");
TClonesArray *branchElectronPF = treeReaderElectron->UseBranch("ElectronPF");
TClonesArray *branchParticleMuon = treeReaderMuon->UseBranch("Particle");
TClonesArray *branchTrackMuon = treeReaderMuon->UseBranch("Track");
TClonesArray *branchMuon = treeReaderMuon->UseBranch("Muon");
TClonesArray *branchParticlePion = treeReaderPion->UseBranch("Particle");
TClonesArray *branchTrackPion = treeReaderPion->UseBranch("Track");
TClonesArray *branchPion = treeReaderPion->UseBranch("Pion");
TClonesArray *branchParticlePhoton = treeReaderPhoton->UseBranch("Particle");
TClonesArray *branchTowerPhoton = treeReaderPhoton->UseBranch("Tower");
TClonesArray *branchPhoton = treeReaderPhoton->UseBranch("Photon");
TClonesArray *branchParticleNeutralHadron = treeReaderNeutralHadron->UseBranch("Particle");
TClonesArray *branchTowerNeutralHadron = treeReaderNeutralHadron->UseBranch("Tower");
TClonesArray *branchGenJet = treeReaderJet->UseBranch("GenJet");
TClonesArray *branchParticleJet = treeReaderJet->UseBranch("Particle");
TClonesArray *branchPFJet = treeReaderJet->UseBranch("PFJet");
TClonesArray *branchCaloJet = treeReaderJet->UseBranch("CaloJet");
TClonesArray *branchJet = treeReaderJet->UseBranch("Jet");
TClonesArray *branchParticleBJet = treeReaderBJet->UseBranch("Particle");
TClonesArray *branchPFBJet = treeReaderBJet->UseBranch("Jet");
TClonesArray *branchParticleCJet = treeReaderCJet->UseBranch("Particle");
TClonesArray *branchPFCJet = treeReaderCJet->UseBranch("Jet");
TClonesArray *branchParticleTauJet = treeReaderTauJet->UseBranch("Particle");
TClonesArray *branchPFTauJet = treeReaderTauJet->UseBranch("Jet");
TClonesArray *branchGenScalarHT = treeReaderJet->UseBranch("GenScalarHT");
TClonesArray *branchMet = treeReaderJet->UseBranch("PFMissingET");
TClonesArray *branchCaloMet = treeReaderJet->UseBranch("CaloMissingET");
std::vector colors;
colors.push_back(kBlack);
colors.push_back(kBlue);
colors.push_back(kRed);
colors.push_back(kGreen+1);
colors.push_back(kMagenta+1);
colors.push_back(kOrange);
std::vector markerStyles;
markerStyles.push_back(20);
markerStyles.push_back(21);
markerStyles.push_back(22);
markerStyles.push_back(23);
markerStyles.push_back(33);
markerStyles.push_back(34);
TString pdfOutput(outputFile);
pdfOutput.ReplaceAll(".root", ".pdf");
TString figPath = inputFilePion;
figPath.ReplaceAll(".root", "");
figPath.ReplaceAll("root", "www/fig");
Int_t lastSlash = figPath.Last('/');
Int_t sizePath = figPath.Length();
figPath.Remove(lastSlash+1,sizePath);
TString header = pdfOutput;
header.ReplaceAll(".pdf", "");
header.ReplaceAll("validation_", "");
lastSlash = header.Last('/');
sizePath = header.Length();
header.Remove(0,lastSlash+1);
TString vrs(version);
TPaveText *pave = new TPaveText(0.0, 0.89, 0.94, 0.94,"NDC");
pave->SetTextAlign(30);
pave->SetTextFont(132);
pave->SetBorderSize(0);
pave->SetShadowColor(0);
pave->SetFillColor(0);
pave->SetFillStyle(0);
pave->AddText("Delphes "+vrs+" - "+header);
TString s_etaMin, s_etaMax, s_eta, s_pt, s_e;
Double_t ptMin = 1.;
Double_t ptMax = 50000.;
Double_t etaMin = -6.;
Double_t etaMax = 6.;
std::vector ptVals;
ptVals.push_back(1.);
ptVals.push_back(10.);
ptVals.push_back(100.);
ptVals.push_back(1000.);
ptVals.push_back(10000.);
std::vector etaVals;
etaVals.push_back(0.);
etaVals.push_back(1.5);
etaVals.push_back(2.5);
etaVals.push_back(4.0);
etaVals.push_back(6.0);
const int n_etabins = etaVals.size()-1;
const int n_ptbins = ptVals.size();
//////////////////////////
// Tracking performance //
//////////////////////////
// --------- Pion Tracks --------- //
TMultiGraph *mg_trkpi_res_pt = new TMultiGraph("","");
TMultiGraph *mg_trkpi_eff_pt = new TMultiGraph("","");
TMultiGraph *mg_trkpi_res_eta = new TMultiGraph("","");
TMultiGraph *mg_trkpi_eff_eta = new TMultiGraph("","");
TLegend *leg_trkpi_res_pt = new TLegend(0.55,0.22,0.90,0.48);
TLegend *leg_trkpi_eff_pt = (TLegend*)leg_trkpi_res_pt->Clone();
TLegend *leg_trkpi_res_eta = (TLegend*)leg_trkpi_res_pt->Clone();
TLegend *leg_trkpi_eff_eta = (TLegend*)leg_trkpi_res_eta->Clone();
TGraphErrors *gr_trkpi_res_pt = new TGraphErrors[n_etabins];
TGraphErrors *gr_trkpi_eff_pt = new TGraphErrors[n_etabins];
TGraphErrors *gr_trkpi_res_eta = new TGraphErrors[n_ptbins];
TGraphErrors *gr_trkpi_eff_eta = new TGraphErrors[n_ptbins];
TH1D* h_trkpi_eff_pt, *h_trkpi_eff_eta;
std::vector *plots_trkpi_res_pt = new std::vector[n_etabins];
std::vector *plots_trkpi_res_eta = new std::vector[n_ptbins];
// loop over eta bins
for (k = 0; k < etaVals.size()-1; k++)
{
HistogramsCollection(&plots_trkpi_res_pt[k], TMath::Log10(ptMin), TMath::Log10(ptMax), "trkpi");
GetPtres