ATLAS Offline Software
Loading...
Searching...
No Matches
DefParamPullPlots.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3*/
4
9
10namespace Trk {
11 void
13 Pull_d0 = nullptr;
14 Pull_z0 = nullptr;
15 Pull_theta = nullptr;
16 Pull_phi = nullptr;
17 Pull_qOverP = nullptr;
18
19 Pull_d0_vs_pt = nullptr;
20 Pull_z0_vs_pt = nullptr;
21 Pull_theta_vs_pt = nullptr;
22 Pull_phi_vs_pt = nullptr;
23 Pull_qOverP_vs_pt = nullptr;
24 }
25
26 void
28 Pull_d0 =
29 Book1D("Pull" + m_sType + "_d0", "d0 Pull " + m_sType + ";d0 (track - truth) / error;Entries", 100, -5, 5);
30 Pull_z0 =
31 Book1D("Pull" + m_sType + "_z0", "z0 Pull " + m_sType + ";z0 (track - truth) / error;Entries", 100, -5, 5);
32 Pull_phi = Book1D("Pull" + m_sType + "_phi", "phi Pull " + m_sType + ";phi (track - truth) / error;Entries", 100,
33 -5, 5);
34 Pull_theta = Book1D("Pull" + m_sType + "_theta", "theta Pull " + m_sType + ";theta (track - truth) / error;Entries",
35 100, -5, 5);
36 Pull_qOverP = Book1D("Pull" + m_sType + "_qOverP",
37 "qOverP Pull " + m_sType + ";Q/P (track - truth) / error;Entries", 100, -5, 5);
38
39 Double_t Pull_ptBins[8] = {
40 5, 10, 25, 40, 60, 100, 200, 1000
41 };
42 std::pair<Int_t, Double_t *> ptBins(7, Pull_ptBins);
43 Pull_d0_vs_pt = Book2D("Pull" + m_sType + "_d0_vs_pt", m_sType + " d0 Pull vs pt;pt;d0 Pull", ptBins.first,
44 ptBins.second, 160, -8, 8);
45 Pull_z0_vs_pt = Book2D("Pull" + m_sType + "_z0_vs_pt", m_sType + " z0 Pull vs pt;pt;z0 Pull", ptBins.first,
46 ptBins.second, 160, -8, 8);
47 Pull_phi_vs_pt = Book2D("Pull" + m_sType + "_phi_vs_pt", m_sType + " phi Pull vs pt;pt;phi Pull", ptBins.first,
48 ptBins.second, 160, -8, 8);
49 Pull_theta_vs_pt = Book2D("Pull" + m_sType + "_theta_vs_pt", m_sType + " theta Pull vs pt;pt;theta Pull",
50 ptBins.first, ptBins.second, 160, -8, 8);
51 Pull_qOverP_vs_pt = Book2D("Pull" + m_sType + "_qOverP_vs_pt", m_sType + " qOverP Pull vs pt;pt;Q/P Pull",
52 ptBins.first, ptBins.second, 160, -8, 8);
53 }
54
55 void
56 DefParamPullPlots::fill(const xAOD::TrackParticle &trkprt, const xAOD::TruthParticle &truthprt, float weight) {
57 const double d0_trk = trkprt.d0();
58 const double z0_trk = trkprt.z0();
59 const double phi_trk = trkprt.phi0();
60 const double theta_trk = trkprt.theta();
61 const double qoverp_trk = trkprt.qOverP();
62
63#if !defined(ROOTCORE) && defined(STANDALONE)
64 const Perigee &trackPer = trkprt.perigeeParameters();
65 const AmgSymMatrix(5) * locError = trackPer.covariance();
66
67 const double err_d0 = Amg::error(*locError, d0);
68 const double err_z0 = Amg::error(*locError, z0);
69 const double err_phi = Amg::error(*locError, phi0);
70 const double err_theta = Amg::error(*locError, theta);
71 const double err_qoverp = Amg::error(*locError, qOverP);
72#else
73 const AmgSymMatrix(5) locError = trkprt.definingParametersCovMatrix();
74
75 const double err_d0 = Amg::error(locError, 0); // d0
76 const double err_z0 = Amg::error(locError, 1);// z0
77 const double err_phi = Amg::error(locError, 2); // phi
78 const double err_theta = Amg::error(locError, 3); // theta
79 const double err_qoverp = Amg::error(locError, 4); // qOverP
80#endif
81
82 const double d0_truth = 0;
83 const double phi_truth = truthprt.phi();
84 const double theta_truth = 2 * atan(exp(-truthprt.eta()));
85 const double qoverp_truth = truthprt.charge() / (truthprt.pt() * cosh(truthprt.eta()));
86 const double pt_truth = truthprt.pt() / 1000.;
87
88 const double d0_pull = (d0_trk - d0_truth) / err_d0;
89 const double phi_pull = (phi_trk - phi_truth) / err_phi;
90 const double theta_pull = (theta_trk - theta_truth) / err_theta;
91 const double qoverp_pull = (qoverp_trk - qoverp_truth) / err_qoverp;
92
93 if (truthprt.hasProdVtx()) {
94 const double z0_truth = truthprt.prodVtx()->z() - trkprt.vz();
95 const double z0_pull = (z0_trk - z0_truth) / err_z0;
96 Pull_z0->Fill(z0_pull,weight);
97 Pull_z0_vs_pt->Fill(pt_truth, z0_pull,weight);
98 }
99
100 Pull_d0->Fill(d0_pull,weight);
101 Pull_phi->Fill(phi_pull,weight);
102 Pull_theta->Fill(theta_pull,weight);
103 Pull_qOverP->Fill(qoverp_pull,weight);
104
105 Pull_d0_vs_pt->Fill(pt_truth, d0_pull,weight);
106 Pull_phi_vs_pt->Fill(pt_truth, phi_pull,weight);
107 Pull_theta_vs_pt->Fill(pt_truth, theta_pull,weight);
108 Pull_qOverP_vs_pt->Fill(pt_truth, qoverp_pull,weight);
109 }
110}
#define AmgSymMatrix(dim)
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
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::TrackParticle &trkprt, const xAOD::TruthParticle &truthprt, float weight=1.0)
float z0() const
Returns the parameter.
float theta() const
Returns the parameter, which has range 0 to .
const Trk::Perigee & perigeeParameters() const
Returns the Trk::MeasuredPerigee track parameters.
const ParametersCovMatrix_t definingParametersCovMatrix() const
Returns the 5x5 symmetric matrix containing the defining parameters covariance matrix.
float vz() const
The z origin for the parameters.
float d0() const
Returns the parameter.
float qOverP() const
Returns the parameter.
float phi0() const
Returns the parameter, which has range to .
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.
double charge() const
Physical charge.
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.
double error(const Amg::MatrixX &mat, int index)
return diagonal error of the matrix caller should ensure the matrix is symmetric and the index is in ...
Ensure that the ATLAS eigen extensions are properly loaded.
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee
@ phi0
Definition ParamDefs.h:65
@ theta
Definition ParamDefs.h:66
@ qOverP
perigee
Definition ParamDefs.h:67
@ d0
Definition ParamDefs.h:63
@ z0
Definition ParamDefs.h:64
TrackParticle_v1 TrackParticle
Reference the current persistent version:
TruthParticle_v1 TruthParticle
Typedef to implementation.