ATLAS Offline Software
Loading...
Searching...
No Matches
InDetPerfPlot_TRTExtension.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
9
10
12
15
16#include <cmath>
17
18using namespace IDPVM;
19
21 InDetPlotBase(pParent, sDir),
22 m_resolutionMethod(IDPVM::ResolutionHelper::iterRMS_convergence),
61
62
63 {
64}
65
66void
68
69 book(m_fracTRTExtensions_vs_eta, "fracTRTExtensions_vs_eta");
70 book(m_fracTRTExtensions_vs_pt, "fracTRTExtensions_vs_pt");
71 book(m_fracTRTExtensions_vs_mu, "fracTRTExtensions_vs_mu");
72 book(m_fracTRTExtensions_vs_nvertices, "fracTRTExtensions_vs_nvertices");
73
74 book(m_fracTRTExtensions_matched_vs_eta, "fracTRTExtensions_matched_vs_eta");
75 book(m_fracTRTExtensions_matched_vs_pt, "fracTRTExtensions_matched_vs_pt");
76
77 book(m_fracFindableTRTExtensions_vs_eta, "fracFindableTRTExtensions_vs_eta");
78 book(m_fracFindableTRTExtensions_vs_pt, "fracFindableTRTExtensions_vs_pt");
79 book(m_fracFindableTRTExtensions_vs_mu, "fracFindableTRTExtensions_vs_mu");
80 book(m_fracFindableTRTExtensions_vs_nvertices, "fracFindableTRTExtensions_vs_nvertices");
81
82 book(m_fracFindableTRTExtensions_matched_vs_eta, "fracFindableTRTExtensions_matched_vs_eta");
83 book(m_fracFindableTRTExtensions_matched_vs_pt, "fracFindableTRTExtensions_matched_vs_pt");
84
85 book(m_chi2ndofTRTExtensions, "chi2ndofTRTExtensions");
86 book(m_chi2ndofNoTRTExtensions, "chi2ndofNoTRTExtensions");
87
88 book(m_ptresTRTExtensions_vs_eta, "resHelper_eta_pt_TRTExtension");
89 book(m_ptresTRTExtensions_vs_pt, "resHelper_pt_pt_TRTExtension");
90 book(m_ptresNoTRTExtensions_vs_eta, "resHelper_eta_pt_NoTRTExtension");
91 book(m_ptresNoTRTExtensions_vs_pt, "resHelper_pt_pt_NoTRTExtension");
92
93 book(m_reswidthTRTExtensions_vs_eta, "resolution_vs_eta_pt_TRTExtension");
94 book(m_resmeanTRTExtensions_vs_eta, "resmean_vs_eta_pt_TRTExtension");
95 book(m_reswidthTRTExtensions_vs_pt, "resolution_vs_pt_pt_TRTExtension");
96 book(m_resmeanTRTExtensions_vs_pt, "resmean_vs_pt_pt_TRTExtension");
97
98 book(m_reswidthNoTRTExtensions_vs_eta, "resolution_vs_eta_pt_NoTRTExtension");
99 book(m_resmeanNoTRTExtensions_vs_eta, "resmean_vs_eta_pt_NoTRTExtension");
100 book(m_reswidthNoTRTExtensions_vs_pt, "resolution_vs_pt_pt_NoTRTExtension");
101 book(m_resmeanNoTRTExtensions_vs_pt, "resmean_vs_pt_pt_NoTRTExtension");
102
103
104 book(m_ptpullTRTExtensions_vs_eta, "pullHelper_eta_pt_TRTExtension");
105 book(m_ptpullTRTExtensions_vs_pt, "pullHelper_pt_pt_TRTExtension");
106 book(m_ptpullNoTRTExtensions_vs_eta, "pullHelper_eta_pt_NoTRTExtension");
107 book(m_ptpullNoTRTExtensions_vs_pt, "pullHelper_pt_pt_NoTRTExtension");
108
109 book(m_pullwidthTRTExtensions_vs_eta, "pullwidth_vs_eta_pt_TRTExtension");
110 book(m_pullmeanTRTExtensions_vs_eta, "pullmean_vs_eta_pt_TRTExtension");
111 book(m_pullwidthTRTExtensions_vs_pt, "pullwidth_vs_pt_pt_TRTExtension");
112 book(m_pullmeanTRTExtensions_vs_pt, "pullmean_vs_pt_pt_TRTExtension");
113
114 book(m_pullwidthNoTRTExtensions_vs_eta, "pullwidth_vs_eta_pt_NoTRTExtension");
115 book(m_pullmeanNoTRTExtensions_vs_eta, "pullmean_vs_eta_pt_NoTRTExtension");
116 book(m_pullwidthNoTRTExtensions_vs_pt, "pullwidth_vs_pt_pt_NoTRTExtension");
117 book(m_pullmeanNoTRTExtensions_vs_pt, "pullmean_vs_pt_pt_NoTRTExtension");
118
119}
120
121void
123
124 double eta = particle.eta();
125 double pt = particle.pt() / Gaudi::Units::GeV;
126 float chi2 = particle.chiSquared();
127 float ndof = particle.numberDoF();
128 float chi2Overndof = ndof > 0 ? chi2 / ndof : 0;
129
130 uint8_t iTrtHits = 0;
131 uint8_t iTrtOutliers = 0;
132 particle.summaryValue(iTrtHits, xAOD::numberOfTRTHits);
133 particle.summaryValue(iTrtOutliers, xAOD::numberOfTRTOutliers);
134
135 bool hasTRTHits = iTrtHits + iTrtOutliers > 0;
136 // failed extensions will have only outlier hits and iTRThits will be 0
137 bool isTRTExtension = iTrtHits > 0;
138
139 fillHisto(m_fracTRTExtensions_vs_eta, eta, isTRTExtension, weight);
140 fillHisto(m_fracTRTExtensions_vs_pt, pt, isTRTExtension, weight);
141 if (hasTRTHits) {
142 fillHisto(m_fracFindableTRTExtensions_vs_eta, eta, isTRTExtension, weight);
143 fillHisto(m_fracFindableTRTExtensions_vs_pt, pt, isTRTExtension, weight);
144 }
145
146 if(isTRTExtension) fillHisto(m_chi2ndofTRTExtensions, chi2Overndof, weight);
147 else { fillHisto(m_chi2ndofNoTRTExtensions, chi2Overndof, weight); }
148
149}
150
151void
152InDetPerfPlot_TRTExtension::fill(const xAOD::TrackParticle& particle, const float mu, const unsigned int nvertices, float weight) {
153
154 uint8_t iTrtHits = 0;
155 uint8_t iTrtOutliers = 0;
156 particle.summaryValue(iTrtHits, xAOD::numberOfTRTHits);
157 particle.summaryValue(iTrtOutliers, xAOD::numberOfTRTOutliers);
158
159 bool hasTRTHits = iTrtHits + iTrtOutliers > 0;
160 // failed extensions will have only outlier hits and iTRThits will be 0
161 bool isTRTExtension = iTrtHits > 0;
162
163 fillHisto(m_fracTRTExtensions_vs_mu, mu, isTRTExtension, weight);
164 fillHisto(m_fracTRTExtensions_vs_nvertices, nvertices, isTRTExtension, weight);
165 if (hasTRTHits) {
166 fillHisto(m_fracFindableTRTExtensions_vs_mu, mu, isTRTExtension, weight);
167 fillHisto(m_fracFindableTRTExtensions_vs_nvertices, nvertices, isTRTExtension, weight);
168 }
169
170}
171
172void
173InDetPerfPlot_TRTExtension::fill(const xAOD::TrackParticle& particle, const xAOD::TruthParticle& truthParticle, float weight) {
174
175 //Fraction of extended for truth matched tracks
176
177 uint8_t iTrtHits = 0;
178 uint8_t iTrtOutliers = 0;
179 particle.summaryValue(iTrtHits, xAOD::numberOfTRTHits);
180 particle.summaryValue(iTrtOutliers, xAOD::numberOfTRTOutliers);
181
182 bool hasTRTHits = iTrtHits + iTrtOutliers > 0;
183 // failed extensions will have only outlier hits and iTRThits will be 0
184 bool isTRTExtension = iTrtHits > 0;
185
186 //Get pT resolution for TRT extensions versus without
187 const float undefinedValue = -9999;
188 const float smallestAllowableTan = 1e-8;
189 const float sinTheta{std::sin(particle.theta())};
190 const bool saneSineValue = (std::abs(sinTheta) > 1e-8);
191 const float inverseSinTheta = saneSineValue ? 1./sinTheta : undefinedValue;
192 float track_qopt = saneSineValue ? particle.qOverP()*inverseSinTheta : undefinedValue;
193 const float qopterr = std::sqrt(particle.definingParametersCovMatrix()(4, 4)) * inverseSinTheta;
194
195 static const SG::ConstAccessor<float> qOverPAcc("qOverP");
196 static const SG::ConstAccessor<float> thetaAcc("theta");
197 const float truth_qop = qOverPAcc.isAvailable(truthParticle) ? qOverPAcc(truthParticle) : undefinedValue;
198 const float truth_theta = thetaAcc.isAvailable(truthParticle) ? thetaAcc(truthParticle) : undefinedValue;
199 float truth_qopt = std::abs(truth_theta) > 0 ? truth_qop * 1/(std::sin(truth_theta)) : undefinedValue;
200
201 float ptres = (track_qopt - truth_qopt) * ( 1 / truth_qopt);
202 float ptpull = qopterr > smallestAllowableTan ? (track_qopt - truth_qopt) / qopterr : undefinedValue;
203 float pt = truthParticle.pt() / Gaudi::Units::GeV;
204 const float tanHalfTheta = std::tan(truth_theta * 0.5);
205 const bool tanThetaIsSane = std::abs(tanHalfTheta) > smallestAllowableTan;
206 float eta = undefinedValue;
207 if (tanThetaIsSane) eta = -std::log(tanHalfTheta);
208
209
210 if(isTRTExtension){
212 fillHisto(m_ptresTRTExtensions_vs_pt, pt, ptres, weight);
213
215 fillHisto(m_ptpullTRTExtensions_vs_pt, pt, ptpull, weight);
216 } else {
218 fillHisto(m_ptresNoTRTExtensions_vs_pt, pt, ptres, weight);
219
221 fillHisto(m_ptpullNoTRTExtensions_vs_pt, pt, ptpull, weight);
222 }
223
224 fillHisto(m_fracTRTExtensions_matched_vs_eta, eta, isTRTExtension, weight);
225 fillHisto(m_fracTRTExtensions_matched_vs_pt, pt, isTRTExtension, weight);
226 if (hasTRTHits) {
228 fillHisto(m_fracFindableTRTExtensions_matched_vs_pt, pt, isTRTExtension, weight);
229 }
230
231}
232
233
234void
Scalar eta() const
pseudorapidity method
Helper class to provide constant type-safe access to aux data.
IDPVM::ResolutionHelper::methods m_resolutionMethod
void fill(const xAOD::TrackParticle &particle, float weight)
IDPVM::ResolutionHelper m_resolutionHelper
InDetPerfPlot_TRTExtension(InDetPlotBase *pParent, const std::string &dirName)
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 isAvailable(const ELT &e) const
Test to see if this variable exists in the store.
virtual double pt() const override final
The transverse momentum ( ) of the particle.
double chi2(TH1 *h0, TH1 *h1)
Class to retrieve associated truth from a track, implementing a cached response.
TrackParticle_v1 TrackParticle
Reference the current persistent version:
TruthParticle_v1 TruthParticle
Typedef to implementation.
@ numberOfTRTHits
number of TRT hits [unit8_t].
@ numberOfTRTOutliers
number of TRT outliers [unit8_t].