ATLAS Offline Software
Loading...
Searching...
No Matches
InDetPerfPlot_Efficiency.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
8#include "logLinearBinning.h"
9#include "GaudiKernel/SystemOfUnits.h" //for Gaudi::Units
10using namespace IDPVM;
11
12InDetPerfPlot_Efficiency::InDetPerfPlot_Efficiency(InDetPlotBase* pParent, const std::string& sDir, bool doTechEff, bool isITk) :
13 InDetPlotBase(pParent, sDir),
14 m_doTechEff(doTechEff),
15 m_isITk(isITk) {}
16
17void
19 // Drop eta bins larger than 2.5
20 if(!m_isITk){
22 [](float eta) { return eta > 2.5; }),
23 m_eta_bins.end());
24 }
25
26 book(m_efficiency_vs_pteta, "efficiency_vs_pteta");
27 book(m_efficiency_vs_ptTruthMu, "efficiency_vs_ptTruthMu");
28 book(m_efficiency_vs_ptActualMu, "efficiency_vs_ptActualMu");
29
30 book(m_efficiency_vs_etaTruthMu, "efficiency_vs_absEtaTruthMu");
31 book(m_efficiency_vs_etaActualMu, "efficiency_vs_absEtaActualMu");
32
33 for(unsigned int i=0; i<m_eta_bins.size()-1; i++){
36
37 std::string bin_low = std::to_string(m_eta_bins[i]);
38 size_t dotPos = bin_low.find('.');
39 bin_low.resize(dotPos+2);
40 bin_low.replace(dotPos, 1, 1, 'p');
41 std::string bin_up = std::to_string(m_eta_bins[i+1]);
42 dotPos = bin_up.find('.');
43 bin_up.resize(dotPos+2);
44 bin_up.replace(dotPos, 1, 1, 'p');
45
47 "efficiency_vs_truthMu_absEta_"+bin_low+"_"+bin_up);
49 "efficiency_vs_actualMu_absEta_"+bin_low+"_"+bin_up);
50 }
51
52 book(m_efficiency_vs_eta, "efficiency_vs_eta");
53 book(m_efficiency_vs_pt, "efficiency_vs_pt");
54 book(m_efficiency_vs_pt_low, "efficiency_vs_pt_low");
55 book(m_efficiency_vs_pt_high, "efficiency_vs_pt_high");
56 book(m_efficiency_vs_lowpt, "efficiency_vs_lowpt");
57 book(m_efficiency_vs_phi, "efficiency_vs_phi");
58 book(m_efficiency_vs_d0, "efficiency_vs_d0");
59 book(m_efficiency_vs_d0_abs, "efficiency_vs_d0_abs");
60 book(m_efficiency_vs_z0, "efficiency_vs_z0");
61 book(m_efficiency_vs_z0_abs, "efficiency_vs_z0_abs");
62 book(m_efficiency_vs_R, "efficiency_vs_R");
63 book(m_efficiency_vs_Z, "efficiency_vs_Z");
64 book(m_efficiency_vs_truthMu, "efficiency_vs_truthMu");
65 book(m_efficiency_vs_actualMu, "efficiency_vs_actualMu");
66
67 if(m_doTechEff){
68 book(m_technical_efficiency_vs_eta, "technical_efficiency_vs_eta");
69 book(m_technical_efficiency_vs_pt, "technical_efficiency_vs_pt");
70 book(m_technical_efficiency_vs_phi, "technical_efficiency_vs_phi");
71 book(m_technical_efficiency_vs_d0, "technical_efficiency_vs_d0");
72 book(m_technical_efficiency_vs_z0, "technical_efficiency_vs_z0");
73 book(m_technical_efficiency_vs_truthMu, "technical_efficiency_vs_truthMu");
74 book(m_technical_efficiency_vs_actualMu, "technical_efficiency_vs_actualMu");
75 book(m_technical_efficiency_vs_R, "technical_efficiency_vs_R");
76 book(m_technical_efficiency_vs_prodR, "technical_efficiency_vs_prodR");
77 book(m_technical_efficiency_vs_prodR_extended, "technical_efficiency_vs_prodR_extended");
78 }
79
80 book(m_extended_efficiency_vs_d0, "extended_efficiency_vs_d0");
81 book(m_extended_efficiency_vs_d0_abs, "extended_efficiency_vs_d0_abs");
82 book(m_extended_efficiency_vs_z0, "extended_efficiency_vs_z0");
83 book(m_extended_efficiency_vs_z0_abs, "extended_efficiency_vs_z0_abs");
84 book(m_efficiency_vs_prodR, "efficiency_vs_prodR");
85 book(m_efficiency_vs_prodR_extended, "efficiency_vs_prodR_extended");
86 book(m_efficiency_vs_prodZ, "efficiency_vs_prodZ");
87 book(m_efficiency_vs_prodZ_extended, "efficiency_vs_prodZ_extended");
88
89 book(m_TrkRec_eta, "TrkRec_eta");
90 book(m_TrkRec_d0, "TrkRec_d0");
91 book(m_TrkRec_prodR, "TrkRec_prodR");
92 book(m_TrkRec_pT, "TrkRec_pT");
93 book(m_TrkRec_truthMu, "TrkRec_truthMu");
94 book(m_TrkRec_actualMu, "TrkRec_actualMu");
95 book(m_TrkRec_eta_d0, "TrkRec_eta_d0");
96 book(m_TrkRec_eta_prodR, "TrkRec_eta_prodR");
97 book(m_TrkRec_eta_pT, "TrkRec_eta_pT");
98
99 book(m_efficiency_vs_pt_log, "efficiency_vs_pt_log");
100 const TH1* h = m_efficiency_vs_pt_log->GetTotalHistogram();
101 int nbins = h->GetNbinsX();
102 std::vector<double> logptbins = IDPVM::logLinearBinning(nbins, h->GetBinLowEdge(1), h->GetBinLowEdge(nbins + 1), false);
103 m_efficiency_vs_pt_log->SetBins(nbins, logptbins.data());
104}
105
106void
107InDetPerfPlot_Efficiency::fill(const xAOD::TruthParticle& truth, const bool isGood, unsigned int truthMu, float actualMu, float weight) {
108 double eta = truth.eta();
109 double pt = truth.pt() / Gaudi::Units::GeV; // convert MeV to GeV
110 double phi = truth.phi();
111
112 fillHisto(m_efficiency_vs_pteta, pt, eta, isGood, weight);
113 fillHisto(m_efficiency_vs_ptTruthMu, pt, truthMu, isGood, weight);
114 fillHisto(m_efficiency_vs_ptActualMu, pt, actualMu, isGood, weight);
115
116 fillHisto(m_efficiency_vs_etaTruthMu, std::abs(eta), truthMu, isGood, weight);
117 fillHisto(m_efficiency_vs_etaActualMu, std::abs(eta), actualMu, isGood, weight);
118
119 const auto pVal = std::lower_bound(m_eta_bins.begin(), m_eta_bins.end(), std::abs(eta));
120 const int bin = std::distance(m_eta_bins.begin(), pVal) - 1;
121 // Protect against overflow
122 // for LRT config in particular which keeps truth particles up to eta=3.0
123 if(bin<static_cast<int>(m_efficiency_vs_truthMu_eta_bin.size())){
124 fillHisto(m_efficiency_vs_truthMu_eta_bin[bin], truthMu, isGood, weight);
125 fillHisto(m_efficiency_vs_actualMu_eta_bin[bin], actualMu, isGood, weight);
126 }
127
128 fillHisto(m_efficiency_vs_eta, eta, isGood, weight);
129 fillHisto(m_efficiency_vs_pt, pt, isGood, weight);
130 fillHisto(m_efficiency_vs_pt_low, pt, isGood, weight);
131 fillHisto(m_efficiency_vs_pt_high, pt, isGood, weight);
132 fillHisto(m_efficiency_vs_phi, phi, isGood, weight);
133 fillHisto(m_efficiency_vs_pt_log, pt, isGood, weight);
134 fillHisto(m_efficiency_vs_lowpt, pt, isGood, weight);
135
136 static const SG::ConstAccessor<float> d0Acc("d0");
137 static const SG::ConstAccessor<float> z0Acc("z0");
138 static const SG::ConstAccessor<float> prodRAcc("prodR");
139 static const SG::ConstAccessor<float> prodZAcc("prodZ");
140 double d0 = d0Acc(truth);
141 double z0 = z0Acc(truth);
142 double R = prodRAcc(truth);
143 double Z = prodZAcc(truth);
144 fillHisto(m_efficiency_vs_d0, d0, isGood, weight);
145 fillHisto(m_efficiency_vs_d0_abs, std::abs(d0), isGood, weight);
146 fillHisto(m_efficiency_vs_z0, z0, isGood, weight);
147 fillHisto(m_efficiency_vs_z0_abs, std::abs(z0), isGood, weight);
148 fillHisto(m_efficiency_vs_R, R, isGood, weight);
149 fillHisto(m_efficiency_vs_Z, Z, isGood, weight);
150
151 fillHisto(m_extended_efficiency_vs_d0, d0, isGood, weight);
152 fillHisto(m_extended_efficiency_vs_d0_abs, std::abs(d0), isGood, weight);
153 fillHisto(m_extended_efficiency_vs_z0, z0, isGood, weight);
154 fillHisto(m_extended_efficiency_vs_z0_abs, std::abs(z0), isGood, weight);
155 fillHisto(m_efficiency_vs_truthMu, truthMu, isGood, weight);
156 fillHisto(m_efficiency_vs_actualMu, actualMu, isGood, weight);
157
158 fillHisto(m_TrkRec_eta, eta, isGood, weight);
159 fillHisto(m_TrkRec_d0, d0, isGood, weight);
160 fillHisto(m_TrkRec_pT, pt, isGood, weight);
161 fillHisto(m_TrkRec_truthMu, truthMu, isGood, weight);
162 fillHisto(m_TrkRec_actualMu, actualMu, isGood, weight);
163
164 fillHisto(m_TrkRec_eta_d0, eta, d0, isGood, weight);
165 fillHisto(m_TrkRec_eta_pT, eta, pt, isGood, weight);
166
167 if (truth.hasProdVtx()) {
168 const xAOD::TruthVertex* vtx = truth.prodVtx();
169 double prod_rad = vtx->perp();
170 double prod_z = vtx->z();
171 fillHisto(m_efficiency_vs_prodR, prod_rad, isGood, weight);
172 fillHisto(m_efficiency_vs_prodR_extended, prod_rad, isGood, weight);
173 fillHisto(m_efficiency_vs_prodZ, prod_z, isGood, weight);
174 fillHisto(m_efficiency_vs_prodZ_extended, prod_z, isGood, weight);
175
176 fillHisto(m_TrkRec_prodR, prod_rad, isGood, weight);
177
178 fillHisto(m_TrkRec_eta_prodR, eta, prod_rad, isGood, weight);
179 }
180}
181
182void
183InDetPerfPlot_Efficiency::fillTechnicalEfficiency(const xAOD::TruthParticle& truth, const bool isGood, unsigned int truthMu, float actualMu, float weight) {
184 double eta = truth.eta();
185 double pt = truth.pt() / Gaudi::Units::GeV; // convert MeV to GeV
186 double phi = truth.phi();
187 fillHisto(m_technical_efficiency_vs_pteta, pt, eta, isGood, weight);
189 fillHisto(m_technical_efficiency_vs_pt, pt, isGood, weight);
191
192 static const SG::ConstAccessor<float> d0Acc("d0");
193 static const SG::ConstAccessor<float> z0Acc("z0");
194 static const SG::ConstAccessor<float> prodRAcc("prodR");
195 double d0 = d0Acc(truth);
196 double z0 = z0Acc(truth);
197 double R = prodRAcc(truth);
198 fillHisto(m_technical_efficiency_vs_d0, d0, isGood, weight);
199 fillHisto(m_technical_efficiency_vs_z0, z0, isGood, weight);
200 fillHisto(m_technical_efficiency_vs_truthMu, truthMu, isGood, weight);
201 fillHisto(m_technical_efficiency_vs_actualMu, actualMu, isGood, weight);
202 fillHisto(m_technical_efficiency_vs_R, R, isGood, weight);
203
204 if (truth.hasProdVtx()) {
205 const xAOD::TruthVertex* vtx = truth.prodVtx();
206 double prod_rad = vtx->perp();
207 fillHisto(m_technical_efficiency_vs_prodR, prod_rad, isGood, weight);
208 fillHisto(m_technical_efficiency_vs_prodR_extended, prod_rad, isGood, weight);
209 }
210}
211
212
213void
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
Helper class to provide constant type-safe access to aux data.
InDetPerfPlot_Efficiency(InDetPlotBase *pParent, const std::string &dirName, bool doTechEff=false, bool isITk=false)
void fill(const xAOD::TruthParticle &truth, const bool isGood, unsigned int truthMu, float actualMu, float weight)
void fillTechnicalEfficiency(const xAOD::TruthParticle &truth, const bool isGood, unsigned int truthMu, float actualMu, float weight)
std::vector< TEfficiency * > m_efficiency_vs_truthMu_eta_bin
TEfficiency * m_technical_efficiency_vs_prodR_extended
std::vector< TEfficiency * > m_efficiency_vs_actualMu_eta_bin
static void fillHisto(TProfile *pTprofile, const float bin, const float weight, const float weight2=1.0)
void book(Htype *&pHisto, const std::string &histoIdentifier, const std::string &nameOverride="", const std::string &folder="default")
Helper method to book histograms using an identifier string.
InDetPlotBase(InDetPlotBase *pParent, const std::string &dirName)
Constructor taking parent node and directory name for plots.
Helper class to provide constant type-safe access to aux data.
bool hasProdVtx() const
Check for a production vertex on this particle.
virtual double pt() const override final
The transverse momentum ( ) of the particle.
const TruthVertex_v1 * prodVtx() const
The production vertex of this particle.
virtual double eta() const override final
The pseudorapidity ( ) of the particle.
virtual double phi() const override final
The azimuthal angle ( ) of the particle.
float z() const
Vertex longitudinal distance along the beam line form the origin.
float perp() const
Vertex transverse distance from the beam line.
Class to retrieve associated truth from a track, implementing a cached response.
std::vector< double > logLinearBinning(const unsigned int nBins, const double absXmin, const double absXmax, const bool symmetriseAroundZero)
DataModel_detail::iterator< DVL > remove_if(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end, Predicate pred)
Specialization of remove_if for DataVector/List.
TruthVertex_v1 TruthVertex
Typedef to implementation.
Definition TruthVertex.h:15
TruthParticle_v1 TruthParticle
Typedef to implementation.