/*
* 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;
//------------------------------------------------------------------------------
double ptrangemin = 10;
double ptrangemax = 10000;
static const int Nbins = 20;
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;
struct resolPlot
{
TH1 *cenResolHist;
TH1 *fwdResolHist;
double ptmin;
double ptmax;
double xmin;
double xmax;
TString obj;
resolPlot();
resolPlot(double ptdown, double ptup, TString object);
void set(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);
}
void resolPlot::set(double ptdown, double ptup, TString object, double xmin, double xmax)
{
ptmin = ptdown;
ptmax = ptup;
obj = object;
cenResolHist = new TH1D(obj+"_delta_pt_"+Form("%4.2f",ptmin)+"_"+Form("%4.2f",ptmax)+"_cen", obj+"_delta_pt_"+Form("%4.2f",ptmin)+"_"+Form("%4.2f",ptmax)+"_cen", 200, xmin, xmax);
fwdResolHist = new TH1D(obj+"_delta_pt_"+Form("%4.2f",ptmin)+"_"+Form("%4.2f",ptmax)+"_fwd", obj+"_delta_pt_"+Form("%4.2f",ptmin)+"_"+Form("%4.2f",ptmax)+"_fwd", 200, 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);
}
}
//------------------------------------------------------------------------------
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
std::pair GetEff(TClonesArray *branchReco, TClonesArray *branchParticle, TString name, int pdgID, 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 *histGenPtcen = new TH1D(name+" gen spectra Pt",name+" gen spectra cen", Nbins, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax));
TH1D *histRecoPtcen = new TH1D(name+" reco spectra Pt",name+" reco spectra cen", Nbins, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax));
TH1D *histGenPtfwd = new TH1D(name+" gen spectra Eta",name+" gen spectra fwd", Nbins, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax));
TH1D *histRecoPtfwd = new TH1D(name+" reco spectra Eta",name+" reco spectra fwd", Nbins, TMath::Log10(ptrangemin), TMath::Log10(ptrangemax));
histGenPtcen->SetDirectory(0);
histRecoPtcen->SetDirectory(0);
histGenPtfwd->SetDirectory(0);
histRecoPtfwd->SetDirectory(0);
BinLogX(histGenPtcen);
BinLogX(histRecoPtcen);
BinLogX(histGenPtfwd);
BinLogX(histRecoPtfwd);
// 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;
if (particle->PID == pdgID && genMomentum.Pt() > ptrangemin && genMomentum.Pt() < ptrangemax )
{
// 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) 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;
}
}
pt = genMomentum.Pt();
eta = genMomentum.Eta();
if (TMath::Abs(eta) < 1.5)
{
histGenPtcen->Fill(pt);
if(deltaR < 0.3) { histRecoPtcen->Fill(pt); }
}
else if (TMath::Abs(eta) < 2.5)
{
histGenPtfwd->Fill(pt);
if(deltaR < 0.3) { histRecoPtfwd->Fill(pt); }
}
}
}
}
std::pair histos;
histRecoPtcen->Divide(histGenPtcen);
histRecoPtfwd->Divide(histGenPtfwd);
histos.first = histRecoPtcen;
histos.second = histRecoPtfwd;
return histos;
}
template
void GetEres(std::vector *histos, TClonesArray *branchReco, TClonesArray *branchParticle, int pdgID, ExRootTreeReader *treeReader)
{
Long64_t allEntries = treeReader->GetEntries();
cout << "** Computing 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.E();
eta = TMath::Abs(bestGenMomentum.Eta());
for (bin = 0; bin < Nbins; bin++)
{
if(pt > histos->at(bin).ptmin && pt < histos->at(bin).ptmax && eta < 2.5)
{
if (eta < 1.5) {histos->at(bin).cenResolHist->Fill(recoMomentum.E()/bestGenMomentum.E());}
else if (eta < 2.5) {histos->at(bin).fwdResolHist->Fill(recoMomentum.E()/bestGenMomentum.E());}
}
}
}
}
}
}
template
void GetPtres(std::vector *histos, TClonesArray *branchReco, TClonesArray *branchParticle, int pdgID, 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 < 2.5)
{
if (eta < 1.5) {histos->at(bin).cenResolHist->Fill(recoMomentum.Pt()/bestGenMomentum.Pt());}
else if (eta < 2.5) {histos->at(bin).fwdResolHist->Fill(recoMomentum.Pt()/bestGenMomentum.Pt());}
}
}
}
}
}
}
void GetJetsEres(std::vector *histos, TClonesArray *branchJet, TClonesArray *branchGenJet, ExRootTreeReader *treeReader)
{
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.25)
{
pt = genJetMomentum.E();
eta = TMath::Abs(genJetMomentum.Eta());
for (bin = 0; bin < Nbins; bin++)
{
if(pt > histos->at(bin).ptmin && pt < histos->at(bin).ptmax && eta < 1.5)
{
histos->at(bin).cenResolHist->Fill(jetMomentum.E()/bestGenJetMomentum.E());
}
else if(pt > histos->at(bin).ptmin && pt < histos->at(bin).ptmax && eta < 2.5)
{
histos->at(bin).fwdResolHist->Fill(jetMomentum.E()/bestGenJetMomentum.E());
}
}
}
}
}
}
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).cenResolHist->Fill(met->P4().Px());
histos->at(bin).fwdResolHist->Fill(met->P4().Py());
}
}
}
}
}
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 central, bool rms = false)
{
Int_t bin;
Int_t count = 0;
TGraphErrors gr = TGraphErrors(Nbins/2);
Double_t sig = 0;
Double_t sigErr = 0;
for (bin = 0; bin < Nbins; bin++)
{
if (central == true && histos->at(bin).cenResolHist->GetEntries() > 100)
{
std::pair sigvalues = GausFit(histos->at(bin).cenResolHist);
if (rms == true)
{
gr.SetPoint(count,(histos->at(bin).ptmin+histos->at(bin).ptmax)/2.0, sigvalues.second);
gr.SetPointError(count,0, sigvalues.second); // to correct
}
else
{
gr.SetPoint(count,(histos->at(bin).ptmin+histos->at(bin).ptmax)/2.0, sigvalues.first);
gr.SetPointError(count,0, sigvalues.second);
}
count++;
}
else if (central == false && histos->at(bin).fwdResolHist->GetEntries() > 10)
{
std::pair sigvalues = GausFit(histos->at(bin).fwdResolHist);
if (rms == true)
{
gr.SetPoint(count,(histos->at(bin).ptmin+histos->at(bin).ptmax)/2.0, sigvalues.second);
gr.SetPointError(count,0, sigvalues.second); // to correct
}
else
{
gr.SetPoint(count,(histos->at(bin).ptmin+histos->at(bin).ptmax)/2.0, sigvalues.first);
gr.SetPointError(count,0, sigvalues.second);
}
count++;
}
}
return gr;
}
//------------------------------------------------------------------------------
// type 1 : object, 2 : track, 3 : tower
void addGraph(TMultiGraph *mg, TGraphErrors *gr, TLegend *leg, int type)
{
gr->SetLineWidth(2);
switch ( type )
{
case 1:
gr->SetLineColor(objColor);
gr->SetLineStyle(objStyle);
std::cout << "Adding " << gr->GetName() << std::endl;
mg->Add(gr);
leg->AddEntry(gr,"Reco","l");
break;
case 2:
gr->SetLineColor(trackColor);
gr->SetLineStyle(trackStyle);
mg->Add(gr);
leg->AddEntry(gr,"Track","l");
break;
case 3:
gr->SetLineColor(towerColor);
gr->SetLineStyle(towerStyle);
mg->Add(gr);
leg->AddEntry(gr,"Tower","l");
break;
case 0:
gr->SetLineColor(objColor);
gr->SetLineStyle(objStyle);
mg->Add(gr);
break;
default:
std::cout << "wrong type, possibles choices are Object, Track and Tower" << std::endl;
break;
}
}
void addHist(TH1D *h, TLegend *leg, int type)
{
h->SetLineWidth(2);
switch ( type )
{
case 1:
h->SetLineColor(objColor);
h->SetLineStyle(objStyle);
leg->AddEntry(h,"Reco","l");
break;
case 2:
h->SetLineColor(trackColor);
h->SetLineStyle(trackStyle);
leg->AddEntry(h,"Track","l");
break;
case 3:
h->SetLineColor(towerColor);
h->SetLineStyle(towerStyle);
leg->AddEntry(h,"Tower","l");
break;
case 0:
h->SetLineColor(objColor);
h->SetLineStyle(objStyle);
break;
default:
std::cout << "wrong type, possibles choices are Object, Track and Tower" << std::endl;
break;
}
}
void DrawAxis(TMultiGraph *mg, TLegend *leg, double max, int type = 0)
{
mg->SetMinimum(0.);
mg->SetMaximum(max);
mg->GetXaxis()->SetLimits(ptrangemin,ptrangemax);
mg->GetYaxis()->SetTitle("resolution");
if (type == 0) mg->GetXaxis()->SetTitle("E [GeV]");
else mg->GetXaxis()->SetTitle("p_{T} [GeV]");
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->GetYaxis()->SetNdivisions(505);
leg->SetBorderSize(0);
leg->SetShadowColor(0);
leg->SetFillColor(0);
leg->SetFillStyle(0);
gStyle->SetOptTitle(0);
gPad->SetLogx();
gPad->SetBottomMargin(0.2);
gPad->SetLeftMargin(0.2);
gPad->Modified();
gPad->Update();
}
void DrawAxis(TH1D *h, TLegend *leg, int type = 0)
{
h->GetYaxis()->SetRangeUser(0,1.0);
if (type == 0) h->GetXaxis()->SetTitle("E [GeV]");
else h->GetXaxis()->SetTitle("p_{T} [GeV]");
h->GetYaxis()->SetTitle("efficiency");
h->GetYaxis()->SetTitleSize(0.07);
h->GetXaxis()->SetTitleSize(0.07);
h->GetYaxis()->SetLabelSize(0.06);
h->GetXaxis()->SetLabelSize(0.06);
h->GetYaxis()->SetLabelOffset(0.03);
h->GetYaxis()->SetTitleOffset(1.3);
h->GetXaxis()->SetTitleOffset(1.4);
h->GetYaxis()->SetNdivisions(505);
leg->SetBorderSize(0);
leg->SetShadowColor(0);
leg->SetFillColor(0);
leg->SetFillStyle(0);
gStyle->SetOptTitle(0);
gStyle->SetOptStat(0);
gPad->SetBottomMargin(0.2);
gPad->SetLeftMargin(0.2);
gPad->Modified();
gPad->Update();
}
void Validation(const char *inputFileElectron, const char *inputFileMuon, const char *inputFilePhoton, const char *inputFileJet, const char *inputFileBJet, const char *inputFileTauJet, const char *outputFile)
{
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 *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 *chainTauJet = new TChain("Delphes");
chainTauJet->Add(inputFileTauJet);
ExRootTreeReader *treeReaderTauJet = new ExRootTreeReader(chainTauJet);
TClonesArray *branchParticleElectron = treeReaderElectron->UseBranch("Particle");
TClonesArray *branchTrackElectron = treeReaderElectron->UseBranch("Track");
TClonesArray *branchTowerElectron = treeReaderElectron->UseBranch("Tower");
TClonesArray *branchElectron = treeReaderElectron->UseBranch("Electron");
TClonesArray *branchParticleMuon = treeReaderMuon->UseBranch("Particle");
TClonesArray *branchTrackMuon = treeReaderMuon->UseBranch("Track");
TClonesArray *branchMuon = treeReaderMuon->UseBranch("Muon");
TClonesArray *branchParticlePhoton = treeReaderPhoton->UseBranch("Particle");
TClonesArray *branchTowerPhoton = treeReaderPhoton->UseBranch("Tower");
TClonesArray *branchPhoton = treeReaderPhoton->UseBranch("Photon");
TClonesArray *branchGenJet = treeReaderJet->UseBranch("GenJet");
TClonesArray *branchPFJet = treeReaderJet->UseBranch("Jet");
TClonesArray *branchCaloJet = treeReaderJet->UseBranch("CaloJet");
TClonesArray *branchParticleBJet = treeReaderBJet->UseBranch("Particle");
TClonesArray *branchPFBJet = treeReaderBJet->UseBranch("Jet");
TClonesArray *branchParticleTauJet = treeReaderTauJet->UseBranch("Particle");
TClonesArray *branchPFTauJet = treeReaderTauJet->UseBranch("Jet");
TClonesArray *branchScalarHT = treeReaderJet->UseBranch("ScalarHT");
TClonesArray *branchMet = treeReaderJet->UseBranch("MissingET");
///////////////
// Electrons //
///////////////
// Reconstruction efficiency
TString elecs = "Electron";
int elID = 11;
std::pair histos_el = GetEff(branchElectron, branchParticleElectron, "Electron", elID, treeReaderElectron);
// tracking reconstruction efficiency
std::pair histos_eltrack = GetEff