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