ATLAS Offline Software
Loading...
Searching...
No Matches
InDetPerfPlot_Resolution.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
9
11#include <map>
12#include <utility>
13#include <cmath>
14#include "logLinearBinning.h"
17
20
21
22namespace{
23 constexpr float undefinedValue{-9999.};
24 constexpr float smallestAllowableSin(1e-8);
25 constexpr float smallestAllowableTan(1e-8);
26 constexpr float smallestAllowableQoverPt(1e-8);
27}
28
29InDetPerfPlot_Resolution::InDetPerfPlot_Resolution(InDetPlotBase* pParent, const std::string& sDir, bool d0Only, bool hasHGTDReco) : InDetPlotBase(pParent, sDir),
30 m_resolutionMethod(IDPVM::ResolutionHelper::iterRMS_convergence),
31 m_primTrk(false),
32 m_secdTrk(false),
33 m_allTrk(false),
34
35 m_trkP{},
36 m_truetrkP{},
37 m_trkErrP{},
38
39 m_resP{},
40 m_pullP{},
41 m_sigP{},
42
43 m_pull{},
44 m_res{},
45 m_sigma{},
46
53
60
67
74
81
88
95
99 {
100 m_d0Only = d0Only;
101 m_hasHGTDReco = hasHGTDReco;
102
103 TString tsDir = (TString) sDir;
104
105 if (tsDir.Contains("Primary")) {
106 m_primTrk = true; // control if sec/prim from init
107 } else if (tsDir.Contains("Secondary")) {
108 m_secdTrk = true;
109 } else {
110 m_allTrk = true;
111 }
112
113
114
115 std::vector<double> ptBins = IDPVM::logLinearBinning(m_nPtBins, m_ptMin, m_ptMax, false);
116 for (int ipt = 0; ipt <= m_nPtBins; ipt++) {
117 m_PtBins[ipt] = (float) ptBins[ipt];
118 }
119
120 std::vector<double> lowPtBins = IDPVM::logLinearBinning(m_nLowPtBins, m_lowPtMin, m_lowPtMax, false);
121 for (int ipt = 0; ipt <= m_nLowPtBins; ipt++) {
122 m_LowPtBins[ipt] = (float) lowPtBins[ipt];
123 }
124
125
126
127}
128
129void
131
132 for (unsigned int iparam = 0; iparam < NPARAMS; iparam++) {
133 //
134 //1D distributions
135 //
136
137 if(iparam == PT) continue;
138 if(m_d0Only && iparam != D0) continue;
139 if(!m_hasHGTDReco && iparam == TIME) continue;
140
141 book(m_pull[iparam], "pull_" + m_paramProp[iparam]);
142 book(m_res[iparam], "res_" + m_paramProp[iparam]);
143
144 book(m_sigma[iparam], "sigma_" + m_paramProp[iparam]); //New One
145 //
146 //2D Distributions to evaluate resolutions vs eta and pT
147 //
148 book(m_resHelpereta[iparam], "resHelper_eta_" + m_paramProp[iparam]);
149 book(m_resHelperpt[iparam], "resHelper_pt_" + m_paramProp[iparam]);
150 book(m_resHelperlowpt[iparam], "resHelper_lowpt_" + m_paramProp[iparam]);
151 book(m_pullHelpereta[iparam], "pullHelper_eta_" + m_paramProp[iparam]);
152 book(m_pullHelperpt[iparam], "pullHelper_pt_" + m_paramProp[iparam]);
153 book(m_pullHelperlowpt[iparam], "pullHelper_lowpt_" + m_paramProp[iparam]);
154 m_resHelperpt[iparam]->GetXaxis()->Set(m_nPtBins,m_PtBins);
155 m_pullHelperpt[iparam]->GetXaxis()->Set(m_nPtBins,m_PtBins);
156 m_resHelperlowpt[iparam]->GetXaxis()->Set(m_nLowPtBins,m_LowPtBins);
157 m_pullHelperlowpt[iparam]->GetXaxis()->Set(m_nLowPtBins,m_LowPtBins);
158 //
159 //1D Histograms for the final resolution and means
160 //
161 book(m_reswidth_vs_eta[iparam], "resolution_vs_eta_" + m_paramProp[iparam]);
162 book(m_resmean_vs_eta[iparam], "resmean_vs_eta_" + m_paramProp[iparam]);
163 book(m_reswidth_vs_pt[iparam], "resolution_vs_pt_" + m_paramProp[iparam]);
164 book(m_resmean_vs_pt[iparam], "resmean_vs_pt_" + m_paramProp[iparam]);
165 book(m_reswidth_vs_lowpt[iparam], "resolution_vs_lowpt_" + m_paramProp[iparam]);
166 book(m_resmean_vs_lowpt[iparam], "resmean_vs_lowpt_" + m_paramProp[iparam]);
167
168 book(m_pullwidth_vs_eta[iparam], "pullwidth_vs_eta_" + m_paramProp[iparam]);
169 book(m_pullmean_vs_eta[iparam], "pullmean_vs_eta_" + m_paramProp[iparam]);
170 book(m_pullwidth_vs_pt[iparam], "pullwidth_vs_pt_" + m_paramProp[iparam]);
171 book(m_pullmean_vs_pt[iparam], "pullmean_vs_pt_" + m_paramProp[iparam]);
172 book(m_pullwidth_vs_lowpt[iparam], "pullwidth_vs_lowpt_" + m_paramProp[iparam]);
173 book(m_pullmean_vs_lowpt[iparam], "pullmean_vs_lowpt_" + m_paramProp[iparam]);
174
175 m_reswidth_vs_pt[iparam]->GetXaxis()->Set(m_nPtBins, m_PtBins);
176 m_resmean_vs_pt[iparam]->GetXaxis()->Set(m_nPtBins, m_PtBins);
177 m_pullwidth_vs_pt[iparam]->GetXaxis()->Set(m_nPtBins, m_PtBins);
178 m_pullmean_vs_pt[iparam]->GetXaxis()->Set(m_nPtBins, m_PtBins);
179 m_reswidth_vs_lowpt[iparam]->GetXaxis()->Set(m_nLowPtBins, m_LowPtBins);
180 m_resmean_vs_lowpt[iparam]->GetXaxis()->Set(m_nLowPtBins, m_LowPtBins);
181 m_pullwidth_vs_lowpt[iparam]->GetXaxis()->Set(m_nLowPtBins, m_LowPtBins);
182 m_pullmean_vs_lowpt[iparam]->GetXaxis()->Set(m_nLowPtBins, m_LowPtBins);
183
184 //
185 //TProfiles of the sqrt(covii)
186 //
187 book(m_sigma_vs_eta[iparam], "sigma_vs_eta_" + m_paramProp[iparam]);
188 book(m_sigma_vs_pt[iparam], "sigma_vs_pt_" + m_paramProp[iparam]);
189 book(m_sigma_vs_lowpt[iparam], "sigma_vs_lowpt_" + m_paramProp[iparam]);
190
191 //NORA
192 //
193 //Projections only done if high level of detail specified
194 //To keep functionalies to produce projections for varied number of pt / eta bins
195 //keep histogram to the old Book1D methods instead of xml definitions
196 //
197 //
198 //Detailed histograms
199 //
200 if(m_iDetailLevel >= 200){
201 book(m_resHelpereta_pos[iparam], "resHelper_eta_" + m_paramProp[iparam], "resHelper_eta_" + m_paramProp[iparam]+"_posQ");
202 book(m_resHelpereta_neg[iparam], "resHelper_eta_" + m_paramProp[iparam], "resHelper_eta_" + m_paramProp[iparam]+"_negQ");
203 book(m_resHelperpt_pos[iparam], "resHelper_pt_" + m_paramProp[iparam], "resHelper_pt_" + m_paramProp[iparam]+"_posQ");
204 book(m_resHelperpt_neg[iparam], "resHelper_pt_" + m_paramProp[iparam], "resHelper_pt_" + m_paramProp[iparam]+"_negQ");
205 book(m_resHelperlowpt_pos[iparam], "resHelper_lowpt_" + m_paramProp[iparam], "resHelper_lowpt_" + m_paramProp[iparam]+"_posQ");
206 book(m_resHelperlowpt_neg[iparam], "resHelper_lowpt_" + m_paramProp[iparam], "resHelper_lowpt_" + m_paramProp[iparam]+"_negQ");
207
208 //Add log binning
209 m_resHelperpt_pos[iparam]->GetXaxis()->Set(m_nPtBins,m_PtBins);
210 m_resHelperpt_neg[iparam]->GetXaxis()->Set(m_nPtBins,m_PtBins);
211 m_resHelperlowpt_pos[iparam]->GetXaxis()->Set(m_nLowPtBins,m_LowPtBins);
212 m_resHelperlowpt_neg[iparam]->GetXaxis()->Set(m_nLowPtBins,m_LowPtBins);
213
214 //Resolution, Resolution Mean, Pull, Pull Mean
215 book(m_reswidth_vs_eta_pos[iparam], "resolution_vs_eta_" + m_paramProp[iparam], "resolution_vs_eta_" + m_paramProp[iparam] + "_posQ");
216 book(m_reswidth_vs_eta_neg[iparam], "resolution_vs_eta_" + m_paramProp[iparam], "resolution_vs_eta_" + m_paramProp[iparam] + "_negQ");
217 book(m_resmean_vs_eta_pos[iparam], "resmean_vs_eta_" + m_paramProp[iparam], "resmean_vs_eta_" + m_paramProp[iparam] + "_posQ");
218 book(m_resmean_vs_eta_neg[iparam], "resmean_vs_eta_" + m_paramProp[iparam], "resmean_vs_eta_" + m_paramProp[iparam] + "_negQ");
219
220 book(m_reswidth_vs_pt_pos[iparam], "resolution_vs_pt_" + m_paramProp[iparam], "resolution_vs_pt_" + m_paramProp[iparam] + "_posQ");
221 book(m_reswidth_vs_pt_neg[iparam], "resolution_vs_pt_" + m_paramProp[iparam], "resolution_vs_pt_" + m_paramProp[iparam] + "_negQ");
222 book(m_resmean_vs_pt_pos[iparam], "resmean_vs_pt_" + m_paramProp[iparam], "resmean_vs_pt_" + m_paramProp[iparam] + "_posQ");
223 book(m_resmean_vs_pt_neg[iparam], "resmean_vs_pt_" + m_paramProp[iparam], "resmean_vs_pt_" + m_paramProp[iparam] + "_negQ");
224
225 book(m_reswidth_vs_lowpt_pos[iparam], "resolution_vs_lowpt_" + m_paramProp[iparam], "resolution_vs_lowpt_" + m_paramProp[iparam] + "_posQ");
226 book(m_reswidth_vs_lowpt_neg[iparam], "resolution_vs_lowpt_" + m_paramProp[iparam], "resolution_vs_lowpt_" + m_paramProp[iparam] + "_negQ");
227 book(m_resmean_vs_lowpt_pos[iparam], "resmean_vs_lowpt_" + m_paramProp[iparam], "resmean_vs_lowpt_" + m_paramProp[iparam] + "_posQ");
228 book(m_resmean_vs_lowpt_neg[iparam], "resmean_vs_lowpt_" + m_paramProp[iparam], "resmean_vs_lowpt_" + m_paramProp[iparam] + "_negQ");
229
230 m_reswidth_vs_pt_pos[iparam]->GetXaxis()->Set(m_nPtBins, m_PtBins);
231 m_resmean_vs_pt_pos[iparam]->GetXaxis()->Set(m_nPtBins, m_PtBins);
232 m_reswidth_vs_pt_neg[iparam]->GetXaxis()->Set(m_nPtBins, m_PtBins);
233 m_resmean_vs_pt_neg[iparam]->GetXaxis()->Set(m_nPtBins, m_PtBins);
234 m_reswidth_vs_lowpt_pos[iparam]->GetXaxis()->Set(m_nLowPtBins, m_LowPtBins);
235 m_resmean_vs_lowpt_pos[iparam]->GetXaxis()->Set(m_nLowPtBins, m_LowPtBins);
236 m_reswidth_vs_lowpt_neg[iparam]->GetXaxis()->Set(m_nLowPtBins, m_LowPtBins);
237 m_resmean_vs_lowpt_neg[iparam]->GetXaxis()->Set(m_nLowPtBins, m_LowPtBins);
238
239 std::string tmpName, tmpTitle;
240
241 int nPtBins = m_pullHelperpt[iparam]->GetNbinsX();
242 int nLowPtBins = m_pullHelperlowpt[iparam]->GetNbinsX();
243 int nEtaBins = m_pullHelpereta[iparam]->GetNbinsX();
244
245 std::shared_ptr<TH1D> refHistEta { m_pullHelpereta[iparam]->ProjectionY("refEta")};
246 std::shared_ptr<TH1D> refHistPt { m_pullHelperpt[iparam]->ProjectionY("refPt")};
247 std::shared_ptr<TH1D> refHistLowPt { m_pullHelperlowpt[iparam]->ProjectionY("refLowPt")};
248
249 //Projections - VERY expensive, hence only at higher detail levels above 200
250 if(m_iDetailLevel > 200){
251 for (int ibins = 0; ibins < nPtBins; ibins++) {
252 tmpName = "pullProjection_pt_" + m_paramProp[iparam] +"_bin_"+ std::to_string(ibins + 1);
253 tmpTitle = tmpName + "; (" + m_paramProp[iparam] + "^{reco}-" + m_paramProp[iparam] +
254 "^{true})/#sigma_{" + m_paramProp[iparam] + "}";
255 m_pullProjections_vs_pt[iparam][ibins] = Book1D(tmpName, refHistPt.get(), tmpTitle , false);
256
257
258 tmpName = "resProjection_pt_" + m_paramProp[iparam] +"_bin_"+ std::to_string(ibins + 1);
259 tmpTitle = tmpName + "; " + m_paramProp[iparam] + "^{reco}-" + m_paramProp[iparam] + "^{true} ";
260 m_resProjections_vs_pt[iparam][ibins] = Book1D(tmpName, refHistPt.get(), tmpTitle , false);
261
262 }
263 for (int ibins = 0; ibins < nLowPtBins; ibins++) {
264 tmpName = "pullProjection_lowpt_" + m_paramProp[iparam] +"_bin_"+ std::to_string(ibins + 1);
265 tmpTitle = tmpName + "; (" + m_paramProp[iparam] + "^{reco}-" + m_paramProp[iparam] +
266 "^{true})/#sigma_{" + m_paramProp[iparam] + "}";
267 m_pullProjections_vs_lowpt[iparam][ibins] = Book1D(tmpName, refHistLowPt.get(), tmpTitle , false);
268
269
270 tmpName = "resProjection_lowpt_" + m_paramProp[iparam] +"_bin_"+ std::to_string(ibins + 1);
271 tmpTitle = tmpName + "; " + m_paramProp[iparam] + "^{reco}-" + m_paramProp[iparam] + "^{true} ";
272 m_resProjections_vs_lowpt[iparam][ibins] = Book1D(tmpName, refHistLowPt.get(), tmpTitle , false);
273
274 }
275 for (int ibins = 0; ibins < nEtaBins; ibins++) {
276 tmpName = "pullProjection_eta_" + m_paramProp[iparam] +"_bin_"+ std::to_string(ibins + 1);
277 tmpTitle = tmpName + "; (" + m_paramProp[iparam] + "^{reco}-" + m_paramProp[iparam] +
278 "^{true})/#sigma_{" + m_paramProp[iparam] + "}";
279 m_pullProjections_vs_eta[iparam][ibins] = Book1D(tmpName, refHistEta.get(), tmpTitle , false);
280
281
282 tmpName = "resProjection_eta_" + m_paramProp[iparam] +"_bin_"+ std::to_string(ibins + 1);
283 tmpTitle = tmpName + "; " + m_paramProp[iparam] + "^{reco}-" + m_paramProp[iparam] + "^{true} ";
284 m_resProjections_vs_eta[iparam][ibins] = Book1D(tmpName, refHistEta.get(), tmpTitle , false);
285 }
286 }
287 }
288 //
289 //End of saving resolution and pull residual binnings
290 //
291
292 }
293
294}
295
296void
297InDetPerfPlot_Resolution::fill(const xAOD::TrackParticle& trkprt, const xAOD::TruthParticle& truthprt, float weight) {
298
299 int isPrimTrk = 0;
300 int isSecdTrk = 0;
301
302 if (HepMC::is_simulation_particle(&truthprt)) {
303 isSecdTrk = 1;
304 } else {
305 if (HepMC::uniqueID(truthprt) > 0) isPrimTrk = 1;
306 }
307
308 // Move on to the next track incase the wrong track category
309 if (!isPrimTrk && !isSecdTrk) {
310 return;
311 }
312 if (!isPrimTrk && m_primTrk) {
313 return;
314 }
315 if (!isSecdTrk && m_secdTrk) {
316 return;
317 }
318
319 // Get track paramters for truth and reco tracks for easy access
320 // store all in track parameter maps trkP/truetrkP
321 getTrackParameters(trkprt);
322 getTrackParameters(truthprt);
324 getPlots(weight);
325
326}
327
328
329void
331 for (int iParams = 0; iParams < NPARAMS; iParams++) m_truetrkP[iParams] = 0.;
332 getTrackParameters(trkprt);
334 getPlots(weight, false);
335}
336
337
338
339void
340InDetPerfPlot_Resolution::getPlots(float weight, bool useTruthKin) {
341 float pt = useTruthKin ? m_truetrkP[PT] : m_trkP[PT];
342 float theta = useTruthKin ? m_truetrkP[THETA] : m_trkP[THETA];
343 const float tanHalfTheta = std::tan(theta * 0.5);
344 const bool tanThetaIsSane = std::abs(tanHalfTheta) > smallestAllowableTan;
345 float eta = undefinedValue;
346 if (tanThetaIsSane) eta = -std::log(tanHalfTheta);
347 for (unsigned int iparam = 0; iparam < NPARAMS; iparam++) {
348 if(iparam == PT) continue;
349 if(m_d0Only && iparam != D0) continue;
350 if(!m_hasHGTDReco && iparam == TIME) continue;
351 m_pull[iparam]->Fill(m_pullP[iparam], weight);
352 m_res[iparam]->Fill(m_resP[iparam], weight);
353 m_sigma[iparam]->Fill(m_sigP[iparam], weight);
354 m_sigma_vs_eta[iparam]->Fill(eta, m_sigP[iparam], weight);
355 m_sigma_vs_pt[iparam]->Fill(pt, m_sigP[iparam], weight);
356 m_sigma_vs_lowpt[iparam]->Fill(pt, m_sigP[iparam], weight);
357 m_resHelpereta[iparam]->Fill(eta, m_resP[iparam], weight);
358 m_resHelperpt[iparam]->Fill(pt, m_resP[iparam], weight);
359 m_resHelperlowpt[iparam]->Fill(pt, m_resP[iparam], weight);
360 m_pullHelperpt[iparam]->Fill(pt, m_pullP[iparam], weight);
361 m_pullHelperlowpt[iparam]->Fill(pt, m_pullP[iparam], weight);
362 m_pullHelpereta[iparam]->Fill(eta, m_pullP[iparam], weight);
363
364 if(m_iDetailLevel >= 200){
365 if (m_trkP[QOVERPT] >= 0.) {
366 m_resHelpereta_pos[iparam]->Fill(eta, m_resP[iparam], weight);
367 m_resHelperpt_pos[iparam]->Fill(pt, m_resP[iparam], weight);
368 m_resHelperlowpt_pos[iparam]->Fill(pt, m_resP[iparam], weight);
369 }
370 if (m_trkP[QOVERPT] < 0.) {
371 m_resHelpereta_neg[iparam]->Fill(eta, m_resP[iparam], weight);
372 m_resHelperpt_neg[iparam]->Fill(pt, m_resP[iparam], weight);
373 m_resHelperlowpt_neg[iparam]->Fill(pt, m_resP[iparam], weight);
374 }
375 }
376
377 }
378
379}
380
381void
383 for (unsigned int iparam = 0; iparam < NPARAMS; iparam++) {
384 m_resP[iparam] = m_trkP[iparam] - m_truetrkP[iparam];
385 if(iparam==PHI) m_resP[iparam] = xAOD::P4Helpers::deltaPhi(m_trkP[iparam], m_truetrkP[iparam]);
386 m_sigP[iparam] = m_trkErrP[iparam];
387 (m_sigP[iparam] != 0) ? m_pullP[iparam] = m_resP[iparam] / m_sigP[iparam] : m_pullP[iparam] = undefinedValue;
388 }
389 bool saneQoverPt = std::abs(m_truetrkP[QOVERPT]) > smallestAllowableQoverPt;
390 if (not saneQoverPt){
391 m_resP[QOVERPT] = undefinedValue;
392 } else {
393 m_resP[QOVERPT] = (m_trkP[QOVERPT] - m_truetrkP[QOVERPT]) * (1 / m_truetrkP[QOVERPT]); // relative q/pt resolution
394 }
395 bool saneQoverPtrec = std::abs(m_trkP[QOVERPT]) > smallestAllowableQoverPt;
396 if (not saneQoverPtrec){
397 m_sigP[QOVERPT] = undefinedValue;
398 } else {
399 m_sigP[QOVERPT] *= 1. / std::abs(m_trkP[QOVERPT]); // relative q/pt error
400 }
401
402 // Missing time information - fill with undefined values
403 if (m_truetrkP[TIME] < -9000. || m_trkP[TIME] < -9000.) {
404 m_resP[TIME] = undefinedValue;
405 m_pullP[TIME] = undefinedValue;
406 }
407 if (m_trkErrP[TIME] < -9000.) {
408 m_sigP[TIME] = undefinedValue;
409 m_pullP[TIME] = undefinedValue;
410 }
411}
412
413void
415
416 m_trkP[D0] = trkprt.d0();
417 m_trkP[Z0] = trkprt.z0();
418 m_trkP[QOVERP] = trkprt.qOverP() * Gaudi::Units::GeV;
419 const float sinTheta{std::sin(trkprt.theta())};
420 const float cosTheta{std::cos(trkprt.theta())};
421 const bool saneSineValue = (std::abs(sinTheta) > smallestAllowableSin);
422 const float inverseSinTheta = saneSineValue ? 1./sinTheta : undefinedValue;
423 m_trkP[QOVERPT] = saneSineValue ?
424 trkprt.qOverP()*inverseSinTheta * Gaudi::Units::GeV : undefinedValue;
425 m_trkP[THETA] = trkprt.theta();
426 m_trkP[PHI] = trkprt.phi0();
427 m_trkP[PT] = trkprt.pt() / Gaudi::Units::GeV;
428 m_trkP[Z0SIN] = trkprt.z0() * sinTheta;
429
430 // Track fit errors
431 m_trkErrP[D0] = trkprt.definingParametersCovMatrix()(0, 0) < 0 ? 0 : std::sqrt(trkprt.definingParametersCovMatrix()(0, 0));
432 m_trkErrP[Z0] = trkprt.definingParametersCovMatrix()(1, 1) < 0 ? 0 : std::sqrt(trkprt.definingParametersCovMatrix()(1, 1));
433 m_trkErrP[PHI] = trkprt.definingParametersCovMatrix()(2, 2) < 0 ? 0 : std::sqrt(trkprt.definingParametersCovMatrix()(2, 2));
434 m_trkErrP[THETA] = trkprt.definingParametersCovMatrix()(3, 3) < 0 ? 0 : std::sqrt(trkprt.definingParametersCovMatrix()(3, 3));
435 m_trkErrP[QOVERP] = trkprt.definingParametersCovMatrix()(4, 4) < 0 ? 0 : std::sqrt(trkprt.definingParametersCovMatrix()(4, 4)) * Gaudi::Units::GeV;
436 float qOverPT_err2 = trkprt.definingParametersCovMatrix()(4, 4) * std::pow(inverseSinTheta, 2)
437 + trkprt.definingParametersCovMatrix()(3, 3) * std::pow(trkprt.qOverP() * cosTheta * std::pow(inverseSinTheta, 2), 2)
438 - 2 * trkprt.qOverP() * cosTheta * trkprt.definingParametersCovMatrix()(3, 4) * std::pow(inverseSinTheta, 3);
439 m_trkErrP[QOVERPT] = qOverPT_err2 > 0 ? std::sqrt(qOverPT_err2) * Gaudi::Units::GeV : 0.;
440 float z0sin2 = std::pow(m_trkErrP[Z0] * std::sin(m_trkP[THETA]), 2) + std::pow(m_trkP[Z0] * m_trkErrP[THETA] * std::cos(m_trkP[THETA]), 2) + 2 * m_trkP[Z0] * std::sin(m_trkP[THETA]) * std::cos(m_trkP[THETA]) * trkprt.definingParametersCovMatrix()(1, 3); // Fixing incorrect formula for z0sin error
441 m_trkErrP[Z0SIN] = z0sin2 < 0 ? 0 : std::sqrt(z0sin2);
442
443
444 // Get error on pT, taken from xAOD::TrackingHelpers.pTErr() but this function only works on a pointer input...
445 if (trkprt.definingParametersCovMatrixVec().size() < 15) {
446 throw std::runtime_error(
447 "TrackParticle without covariance matrix for defining parameters or the covariance matrix is wrong dimensionality.");
448 }
449 if (std::abs(trkprt.qOverP()) < 0) {
450 m_trkErrP[PT] = 0.0;
451 throw std::runtime_error("q/p is zero");
452 } else {
453 double pt = trkprt.pt();
454 double diff_qp = -pt / std::abs(trkprt.qOverP());
455 double diff_theta = trkprt.theta() == 0 ? 0 : pt / std::tan(trkprt.theta());
456 const std::vector<float>& cov = trkprt.definingParametersCovMatrixVec();
457 double pt_err2 = diff_qp * (diff_qp * cov[14] + diff_theta * cov[13]) + diff_theta * diff_theta * cov[9];
458 m_trkErrP[PT] = pt_err2 < 0 ? 0 : std::sqrt(pt_err2) / Gaudi::Units::GeV;
459 }
460
461 // Initialize track time and uncertainty to undefined
462 m_trkP[TIME] = undefinedValue;
463 m_trkErrP[TIME] = undefinedValue;
464 if (m_hasHGTDReco) {
465 static const SG::Accessor<uint8_t> accValidTime("hasValidTime");
466 static const SG::Accessor<float> accTime("time");
467 static const SG::Accessor<float> accTimeRes("timeResolution");
468 if(accValidTime.isAvailable(trkprt) && accTime.isAvailable(trkprt) && accTimeRes.isAvailable(trkprt)) {
469 if (trkprt.hasValidTime()) {
470 m_trkP[TIME] = trkprt.time();
471 m_trkErrP[TIME] = trkprt.timeResolution();
472 }
473 }
474 }
475
476}
477
478
479 void
481 // "d0", "z0", "phi", "theta" , "qOverP"
482 // Perigee for truth particles are in aux container
483 for (int iParams = 0; iParams < NPARAMS; iParams++) {
484 m_truetrkP[iParams] = undefinedValue;
486 if (acc.isAvailable(truthprt)) {
487 m_truetrkP[iParams] = acc(truthprt);
488 }
489 }
490 //rescale Pt
491 m_truetrkP[PT] = truthprt.pt() / Gaudi::Units::GeV;
492 //special cases
493 //need both theta and qOverP for qOverPT
494 //we didnt check qOverP yet, since the closest named variable (strangely; see header) is ptqopt
495 static const SG::ConstAccessor<float> qOverPAcc("qOverP");
496 const float qOverP = qOverPAcc.isAvailable(truthprt) ? qOverPAcc(truthprt) : undefinedValue;
497 if ((qOverP != undefinedValue) and (m_truetrkP[THETA] != undefinedValue)){
498 const float sinTheta =std::sin(m_truetrkP[THETA]);
499 if (std::abs(sinTheta) > smallestAllowableSin){
500 m_truetrkP[QOVERPT] = qOverP / sinTheta * Gaudi::Units::GeV;
501 } else {
502 m_truetrkP[QOVERPT] = undefinedValue;
503 }
504 }
505 m_truetrkP[QOVERP] = qOverP * Gaudi::Units::GeV;
506 //need both z0 and theta for z0sinTheta
507 if ((m_truetrkP[Z0] != undefinedValue) and (m_truetrkP[THETA] != undefinedValue)){
508 const float sinTheta =std::sin(m_truetrkP[THETA]);
509 m_truetrkP[Z0SIN] = m_truetrkP[Z0] * sinTheta;
510 } else {
511 m_truetrkP[Z0SIN] = undefinedValue;
512 }
513
514 // Even if we don't have a track time, we always have a truth particle time
515 static const SG::Accessor<float> accTruthTime("time");
516 m_truetrkP[TIME] = accTruthTime.isAvailable(truthprt) ? accTruthTime(truthprt) : undefinedValue;
517
518}
519
520
521void
523
524 bool saveProjections = (m_iDetailLevel > 200);
525 for (unsigned int iparam = 0; iparam < NPARAMS; iparam++) {
526 if(iparam == PT) continue;
527 if(m_d0Only && iparam != D0) continue;
528 if(!m_hasHGTDReco && iparam == TIME) continue;
529 //
530 //Only save vert detailed information if... high detail level
531 //Reduces output for ART / PhysVal
532 //
533 if(m_iDetailLevel >= 200){
534 m_resolutionHelper.makeResolutions(m_resHelpereta[iparam], m_reswidth_vs_eta[iparam], m_resmean_vs_eta[iparam],
535 (saveProjections ? m_resProjections_vs_eta[iparam] : nullptr), saveProjections, m_resolutionMethod);
536 m_resolutionHelper.makeResolutions(m_resHelperpt[iparam], m_reswidth_vs_pt[iparam], m_resmean_vs_pt[iparam],
537 (saveProjections ? m_resProjections_vs_pt[iparam] : nullptr), saveProjections, m_resolutionMethod);
538 m_resolutionHelper.makeResolutions(m_resHelperlowpt[iparam], m_reswidth_vs_lowpt[iparam], m_resmean_vs_lowpt[iparam],
539 (saveProjections ? m_resProjections_vs_lowpt[iparam] : nullptr), saveProjections, m_resolutionMethod);
540 m_resolutionHelper.makeResolutions(m_pullHelperpt[iparam], m_pullwidth_vs_pt[iparam], m_pullmean_vs_pt[iparam],
541 (saveProjections ? m_pullProjections_vs_pt[iparam] : nullptr), saveProjections, m_resolutionMethod);
542 m_resolutionHelper.makeResolutions(m_pullHelperlowpt[iparam], m_pullwidth_vs_lowpt[iparam], m_pullmean_vs_lowpt[iparam],
543 (saveProjections ? m_pullProjections_vs_lowpt[iparam] : nullptr), saveProjections, m_resolutionMethod);
544 m_resolutionHelper.makeResolutions(m_pullHelpereta[iparam], m_pullwidth_vs_eta[iparam], m_pullmean_vs_eta[iparam],
545 (saveProjections ? m_pullProjections_vs_eta[iparam] : nullptr), saveProjections, m_resolutionMethod);
552 } else {
554 m_resolutionHelper.makeResolutions(m_resHelperpt[iparam], m_reswidth_vs_pt[iparam], m_resmean_vs_pt[iparam], m_resolutionMethod);
559 }
560
561 }
562
563}
Scalar eta() const
pseudorapidity method
Scalar theta() const
theta method
Helper class to provide constant type-safe access to aux data.
IDPVM::ResolutionHelper::methods m_resolutionMethod
TH1 * m_resProjections_vs_eta[NPARAMS][m_nEtaBins]
TH1 * m_pullProjections_vs_lowpt[NPARAMS][m_nLowPtBins]
TH1 * m_resProjections_vs_pt[NPARAMS][m_nPtBins]
void fill(const xAOD::TrackParticle &trkprt, const xAOD::TruthParticle &truthprt, float weight)
TH1 * m_pullProjections_vs_eta[NPARAMS][m_nEtaBins]
void getTrackParameters(const xAOD::TruthParticle &truthprt)
void getPlots(float weight=1.0, bool useTruthKin=true)
float m_LowPtBins[m_nLowPtBins+1]
TH1 * m_pullProjections_vs_pt[NPARAMS][m_nPtBins]
IDPVM::ResolutionHelper m_resolutionHelper
InDetPerfPlot_Resolution(InDetPlotBase *pParent, const std::string &dirName, bool d0Only=false, bool hasHGTDReco=false)
TH1 * m_resProjections_vs_lowpt[NPARAMS][m_nLowPtBins]
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.
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
int m_iDetailLevel
Definition PlotBase.h:101
Helper class to provide type-safe access to aux data.
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.
float z0() const
Returns the parameter.
float timeResolution() const
Returns the time resolution.
float theta() const
Returns the parameter, which has range 0 to .
float time() const
Returns the time.
const ParametersCovMatrix_t definingParametersCovMatrix() const
Returns the 5x5 symmetric matrix containing the defining parameters covariance matrix.
uint8_t hasValidTime() const
Returns whether or not the track has a valid time.
float d0() const
Returns the parameter.
float qOverP() const
Returns the parameter.
std::vector< float > definingParametersCovMatrixVec() const
Returns the length 6 vector containing the elements of defining parameters covariance matrix.
virtual double pt() const override final
The transverse momentum ( ) of the particle.
float phi0() const
Returns the parameter, which has range to .
virtual double pt() const override final
The transverse momentum ( ) of the particle.
int uniqueID(const T &p)
bool is_simulation_particle(const T &p)
Method to establish if a particle (or barcode) was created during the simulation (TODO update to be s...
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)
double deltaPhi(double phiA, double phiB)
delta Phi in range [-pi,pi[
TrackParticle_v1 TrackParticle
Reference the current persistent version:
TruthParticle_v1 TruthParticle
Typedef to implementation.