ATLAS Offline Software
Loading...
Searching...
No Matches
ResoTriggerMuonPlots.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
7
8#include <utility>
9
10using namespace xAOD::P4Helpers;
11
12
13ResoTriggerMuonPlots::ResoTriggerMuonPlots(PlotBase* pParent, const std::string& sDir, std::string sType) :
14 PlotBase(pParent, sDir),
15 m_sType(std::move(sType)),
16 m_pt_slices {"0", "25", "55", "100"},
17 m_etaBins {-3., -2.5, -2.4, -1.918, -1.623, -1.348, -1.2329, -1.1479, -1.05, -0.908, -0.791,
18 -0.652, -0.476, -0.324, -0.132, 0, 0.132, 0.324, 0.476, 0.652, 0.791, 0.908,
19 1.05, 1.1479, 1.2329, 1.348, 1.623, 1.918, 2.4, 2.5, 3.}
20{
21 Res_pT = Book1D("Res" + m_sType + "_pT", "Res" + m_sType + "_pT;(1/pT-1/RECOpT)/(1/RECOpT);Entries", 200, -0.25, 0.25);
22 Res_eta = Book1D("Res" + m_sType + "_eta", "Res" + m_sType + "_eta;eta-RECOeta;Entries", 200, -0.02, 0.02);
23 Res_phi = Book1D("Res" + m_sType + "_phi", "Res" + m_sType + "_phi;phi-RECOphi;Entries", 200, -0.005, 0.005);
24
25 Res_pT_vs_pT = Book2D("Res" + m_sType + "_pT_vs_pT", "Res" + m_sType + "_pT vs pT;pT [GeV];(1/pT-1/RECOpT)/(1/RECOpT)", 150, 0.,
26 150., 100, -0.25, 0.25);
28 Book2D("Res" + m_sType + "_eta_vs_pT", "Res" + m_sType + "_eta vs pT;pT [GeV];eta-RECOeta", 150, 0., 150., 200, -0.02, 0.02);
30 Book2D("Res" + m_sType + "_phi_vs_pT", "Res" + m_sType + "_phi vs pT;pT [GeV];phi-RECOphi", 150, 0., 150., 200, -0.005, 0.005);
31 Res_pT_vs_eta.clear();
32 Res_pT_vs_phi.clear();
33 Res_phi_vs_phi.clear();
34 Res_eta_vs_eta.clear();
35
36 for (unsigned int i = 0; i < m_pt_slices.size() - 1; i++) {
37 Res_pT_vs_eta.push_back(Book2D("Res" + m_sType + "_pT_vs_eta_from" + m_pt_slices[i] + "_to" + m_pt_slices[i + 1] + "GeV",
38 "Res" + m_sType + "_pT vs eta from " + m_pt_slices[i] + " GeV to " + m_pt_slices[i + 1] +
39 " GeV; eta;(1/pT-1/RECOpT)/(1/RECOpT)",
40 m_etaBins.size() - 1, &m_etaBins[0], 100, -0.25, 0.25));
41 Res_pT_vs_phi.push_back(Book2D("Res" + m_sType + "_pT_vs_phi_from" + m_pt_slices[i] + "_to" + m_pt_slices[i + 1] + "GeV",
42 "Res" + m_sType + "_pT vs phi from " + m_pt_slices[i] + " GeV to " + m_pt_slices[i + 1] +
43 " GeV;phi;(1/pT-1/RECOpT)/(1/RECOpT)",
44 96, -M_PI, M_PI, 100, -0.25, 0.25));
45 Res_phi_vs_phi.push_back(Book2D("Res" + m_sType + "_phi_vs_phi_from" + m_pt_slices[i] + "_to" + m_pt_slices[i + 1] + "GeV",
46 "Res" + m_sType + "_phi vs phi from " + m_pt_slices[i] + " GeV to " + m_pt_slices[i + 1] +
47 " GeV;phi;(1/pT-1/RECOpT)/(1/RECOpT)",
48 96, -M_PI, M_PI, 100, -0.25, 0.25));
49 Res_eta_vs_eta.push_back(
50 Book2D("Res" + m_sType + "_eta_vs_eta_from" + m_pt_slices[i] + "_to" + m_pt_slices[i + 1] + "GeV",
51 "Res" + m_sType + "_eta vs eta from " + m_pt_slices[i] + " GeV to " + m_pt_slices[i + 1] + " GeV;eta;eta-RECOeta",
52 m_etaBins.size() - 1, &m_etaBins[0], 100, -0.25, 0.25));
53 }
54}
55
56void ResoTriggerMuonPlots::fill(const xAOD::Muon& Trigmu, const xAOD::Muon& Recomu) {
57 // fill plots Res_pt,eta,phi
58 float pt = 0.001 * Recomu.pt();
59 float eta = Recomu.eta();
60 float phi = Recomu.phi();
61 // float respt = (Trigmu.pt() - Recomu.pt())/Recomu.pt();
62 float respt = (1. / Trigmu.pt() - 1. / Recomu.pt()) / (1. / Recomu.pt());
63
64 const float d_phi = deltaPhi(Trigmu.phi(), phi);
65 const float d_eta = Trigmu.eta() - eta;
66 Res_pT->Fill(respt);
67 Res_eta->Fill(d_eta);
68 Res_phi->Fill(d_phi);
69
70 Res_pT_vs_pT->Fill(pt, respt);
71 Res_eta_vs_pT->Fill(pt, d_eta);
72 Res_phi_vs_pT->Fill(pt, d_phi);
73
74 for (unsigned int i = 0; i < m_pt_slices.size() - 1; i++) {
75 float pt_min = atof(m_pt_slices[i].c_str());
76 float pt_max = atof(m_pt_slices[i + 1].c_str());
77 if (pt >= pt_min && pt < pt_max) {
78 Res_pT_vs_eta[i]->Fill(eta, respt);
79 Res_pT_vs_phi[i]->Fill(phi, respt);
80 Res_phi_vs_phi[i]->Fill(phi, d_phi);
81 Res_eta_vs_eta[i]->Fill(eta, d_eta);
82 }
83 }
84}
85
87 float pt = 0.001 * Recomu.pt();
88 float eta = Recomu.eta();
89 float phi = Recomu.phi();
90 float respt = 0;
91 if (L2SAmu.pt() > 0.) respt = (1. / (1000. * L2SAmu.pt()) - 1. / Recomu.pt()) / (1. / Recomu.pt());
92 if (L2SAmu.pt() < 0.) respt = (1. / (-1000. * L2SAmu.pt()) - 1. / Recomu.pt()) / (1. / Recomu.pt());
93
94 const float d_eta = L2SAmu.eta() - eta;
95 const float d_phi = deltaPhi(L2SAmu.phi(), phi);
96 Res_pT->Fill(respt);
97 Res_eta->Fill(d_eta);
98 Res_phi->Fill(d_phi);
99
100 Res_pT_vs_pT->Fill(pt, respt);
101 Res_eta_vs_pT->Fill(pt, d_eta);
102 Res_phi_vs_pT->Fill(pt, d_phi);
103
104 for (unsigned int i = 0; i < m_pt_slices.size() - 1; i++) {
105 float pt_min = atof(m_pt_slices[i].c_str());
106 float pt_max = atof(m_pt_slices[i + 1].c_str());
107 if (pt >= pt_min && pt < pt_max) {
108 Res_phi_vs_phi[i]->Fill(phi, d_phi);
109 Res_eta_vs_eta[i]->Fill(eta, d_eta);
110 Res_pT_vs_eta[i]->Fill(eta, respt);
111 Res_pT_vs_phi[i]->Fill(phi, respt);
112 }
113 }
114}
115
117 float pt = 0.001 * Recomu.pt();
118 float eta = Recomu.eta();
119 float phi = Recomu.phi();
120 float respt = (1. / (L2CBmu.pt()) - 1. / Recomu.pt()) / (1. / Recomu.pt());
121 ;
122
123 const float d_phi = deltaPhi(L2CBmu.phi() , phi);
124 const float d_eta = L2CBmu.eta() - eta;
125 Res_pT->Fill(respt);
126 Res_eta->Fill(d_eta);
127 Res_phi->Fill(d_phi);
128
129 Res_pT_vs_pT->Fill(pt, respt);
130 Res_eta_vs_pT->Fill(pt, d_eta);
131 Res_phi_vs_pT->Fill(pt, d_phi);
132
133 for (unsigned int i = 0; i < m_pt_slices.size() - 1; i++) {
134 float pt_min = atof(m_pt_slices[i].c_str());
135 float pt_max = atof(m_pt_slices[i + 1].c_str());
136 if (pt >= pt_min && pt < pt_max) {
137 Res_phi_vs_phi[i]->Fill(phi, d_phi);
138 Res_eta_vs_eta[i]->Fill(eta, d_eta);
139 Res_pT_vs_eta[i]->Fill(eta, respt);
140 Res_pT_vs_phi[i]->Fill(phi, respt);
141 }
142 }
143}
#define M_PI
Scalar eta() const
pseudorapidity method
Scalar deltaPhi(const MatrixBase< Derived > &vec) const
Scalar phi() const
phi method
TH1D * Book1D(const std::string &name, const std::string &labels, int nBins, float start, float end, bool prependDir=true)
Book a TH1D histogram.
Definition PlotBase.cxx:94
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
void fill(const xAOD::Muon &Trigmu, const xAOD::Muon &Recomu)
std::vector< double > m_etaBins
std::vector< TH2 * > Res_eta_vs_eta
std::vector< TH2 * > Res_pT_vs_eta
std::vector< TH2 * > Res_pT_vs_phi
std::vector< TH2 * > Res_phi_vs_phi
std::vector< std::string > m_pt_slices
ResoTriggerMuonPlots(PlotBase *pParent, const std::string &sDir, std::string sType="")
virtual double eta() const
The pseudorapidity ( ) of the particle.
virtual double phi() const
The azimuthal angle ( ) of the particle.
virtual double pt() const
The transverse momentum ( ) of the particle.
virtual double pt() const
The transverse momentum ( ) of the particle.
virtual double phi() const
The azimuthal angle ( ) of the particle.
virtual double eta() const
The pseudorapidity ( ) of the particle.
virtual double eta() const
The pseudorapidity ( ) of the particle.
virtual double phi() const
The azimuthal angle ( ) of the particle.
virtual double pt() const
The transverse momentum ( ) of the particle.
STL namespace.
setRcore setEtHad setFside pt
L2CombinedMuon_v1 L2CombinedMuon
Define the latest version of the muon CB class.
Muon_v1 Muon
Reference the current persistent version:
L2StandAloneMuon_v2 L2StandAloneMuon
Define the latest version of the muon SA class.