ATLAS Offline Software
Loading...
Searching...
No Matches
MuonResolutionPlots.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2017 CERN for the benefit of the ATLAS collaboration
3*/
4
8using namespace xAOD::P4Helpers;
9namespace Muon{
10
11MuonResolutionPlots::MuonResolutionPlots(PlotBase* pParent, const std::string& sDir, const std::string& sType, bool doBinnedResolutionPlots):
12 PlotBase(pParent, sDir),
13 m_oResolutionPlots(this, "", sType),
14 m_sType(sType),
15 m_doBinnedResolutionPlots(doBinnedResolutionPlots)
16{
17 const int sizeOf_ptBins=7;
18 Double_t ptBins[sizeOf_ptBins] = {5, 15, 30, 45, 60, 100, 200};
19 const int sizeOf_etaBins=10;
20 Double_t etaBins[sizeOf_etaBins] = {-2.7, -2.5, -1.7, -1.05, -0.1, 0.1, 1.05, 1.7, 2.5, 2.7};
21 const int sizeOf_phiBins=9;
22 Double_t phiBins[sizeOf_phiBins] = {-3.15, -2.4, -1.6, -0.8, 0, 0.8, 1.6, 2.4, 3.15};
23
24 //book 2D hists
26 Res_pT_vs_lowpT = Book2D("Res"+m_sType+"_pT_vs_lowpT","Res"+m_sType+"_pT vs low pT;pT [GeV];(pT-pTtruth)/pTtruth", 16,2.,10., 100,-0.5,0.5);
27 Res_pT_vs_highpT = Book2D("Res"+m_sType+"_pT_vs_highpT","Res"+m_sType+"_pT vs high pT;pT [GeV];(pT-pTtruth)/pTtruth",10,100.,1000., 100,-0.5,0.5);
28
29 Res_pT_vs_pT = Book2D("Res"+m_sType+"_pT_vs_pT","Res"+m_sType+"_pT vs pT;pT [GeV];(pT-pTtruth)/pTtruth",
30 sizeOf_ptBins-1, &ptBins[0], 100,-0.5,0.5);
31 Res_pT_vs_eta = Book2D("Res"+m_sType+"_pT_vs_eta","Res"+m_sType+"_pT vs eta;eta,eta-etatruth",
32 sizeOf_etaBins-1, &etaBins[0], 100,-0.5,0.5);
33 Res_pT_vs_phi = Book2D("Res"+m_sType+"_pT_vs_phi","Res"+m_sType+"_pT vs phi;phi,phi-phitruth",
34 sizeOf_phiBins-1, &phiBins[0], 100,-0.5,0.5);
35
36 Res_eta_vs_pT = Book2D("Res"+m_sType+"_eta_vs_pT","Res"+m_sType+"_eta vs pT;pT [GeV];eta-etatruth",
37 sizeOf_ptBins-1, &ptBins[0], 100,-0.5,0.5);
38 Res_phi_vs_pT = Book2D("Res"+m_sType+"_phi_vs_pT","Res"+m_sType+"_phi vs pT;pT [GeV];phi-phitruth",
39 sizeOf_ptBins-1, &ptBins[0], 100,-0.5,0.5);
40 }
41
42}
43
44 void MuonResolutionPlots::fill(const xAOD::TrackParticle& muPrimaryTrk, const xAOD::TruthParticle& truthMu, float weight) {
45 //fill plots Res_pt,eta,phi
46 m_oResolutionPlots.fill(muPrimaryTrk, truthMu, weight);
47
48 //fill plots Res_pt_vs_pt,eta,phi
49 float pt = 0.001*truthMu.pt();
50 float eta = truthMu.eta();
51 float phi = truthMu.phi();
52 float respt = (muPrimaryTrk.pt() - truthMu.pt())/truthMu.pt();
53
55 if ( pt < Res_pT_vs_lowpT->GetXaxis()->GetXmax() ) Res_pT_vs_lowpT->Fill(pt,respt,weight);
56 if ( pt > Res_pT_vs_highpT->GetXaxis()->GetXmin() ) Res_pT_vs_highpT->Fill(pt,respt,weight);
57
58 Res_pT_vs_pT->Fill(pt,respt,weight);
59 Res_pT_vs_eta->Fill(eta,respt,weight);
60 Res_pT_vs_phi->Fill(phi,respt,weight);
61
62 Res_eta_vs_pT->Fill(pt, muPrimaryTrk.eta()-eta,weight);
63 Res_phi_vs_pT->Fill(pt, deltaPhi(muPrimaryTrk.phi(),phi), weight);
64 }
65}
66
67} // closing namespace Muon
Scalar eta() const
pseudorapidity method
Scalar deltaPhi(const MatrixBase< Derived > &vec) const
Scalar phi() const
phi method
MuonResolutionPlots(PlotBase *pParent, const std::string &sDir, const std::string &sType="", bool doBinnedResolutions=false)
Trk::ResolutionPlots m_oResolutionPlots
void fill(const xAOD::TrackParticle &muontp, const xAOD::TruthParticle &truthprt, float weight=1.0)
PlotBase(PlotBase *parent, const std::string &sDir)
Definition PlotBase.cxx:29
TH2F * Book2D(const std::string &name, const std::string &labels, int nBinsX, float startX, float endX, int nBinsY, float startY, float endY, bool prependDir=true)
Book a TH2F histogram.
Definition PlotBase.cxx:123
virtual double phi() const override final
The azimuthal angle ( ) of the particle (has range to .)
virtual double pt() const override final
The transverse momentum ( ) of the particle.
virtual double eta() const override final
The pseudorapidity ( ) of the particle.
virtual double pt() const override final
The transverse momentum ( ) of the particle.
virtual double eta() const override final
The pseudorapidity ( ) of the particle.
virtual double phi() const override final
The azimuthal angle ( ) of the particle.
NRpcCablingAlg reads raw condition data and writes derived condition data to the condition store.
setRcore setEtHad setFside pt
TrackParticle_v1 TrackParticle
Reference the current persistent version:
TruthParticle_v1 TruthParticle
Typedef to implementation.