ATLAS Offline Software
Loading...
Searching...
No Matches
InDetPerfPlot_TrackParameters.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
10
12
15#include <cmath>
16using namespace IDPVM;
17
19 InDetPlotBase(pParent, sDir), m_isITk{isITk}{
20 //nop
21 //variable initialised at declaration
22}
23
24void
26
27 book(m_reco_d0, "reco_d0");
28 book(m_reco_z0, "reco_z0");
29 book(m_reco_z0sin, "reco_z0sin");
30 book(m_reco_phi, "reco_phi");
31 book(m_reco_theta, "reco_theta");
32 book(m_reco_eta, "reco_eta");
33 book(m_reco_qoverp, "reco_qoverp");
34 book(m_reco_pt, "reco_pt");
35 book(m_reco_lowpt, "reco_lowpt");
36 book(m_reco_chi2, "reco_chi2");
37 book(m_reco_ndof, "reco_ndof");
38 book(m_reco_chi2Overndof, "reco_chi2Overndof");
39 book(m_reco_author, "reco_author");
40 if(m_isITk){
41 book(m_reco_time, "reco_time");
42 book(m_reco_hasValidTime_eff_vs_eta, "reco_hasValidTime_eff_vs_eta");
43 }
44
45 book(m_truth_d0, "truth_d0");
46 book(m_truth_z0, "truth_z0");
47 book(m_truth_z0sin, "truth_z0sin");
48 book(m_truth_phi, "truth_phi");
49 book(m_truth_theta, "truth_theta");
50 book(m_truth_eta, "truth_eta");
51 book(m_truth_qoverp, "truth_qoverp");
52 book(m_truth_pt, "truth_pt");
53 book(m_truth_lowpt, "truth_lowpt");
54 book(m_truth_prodR, "truth_prodR");
55 book(m_truth_prodZ, "truth_prodZ");
56
57 book(m_reco_pt_vs_eta, "reco_pt_vs_eta");
58 book(m_reco_phi_vs_eta, "reco_phi_vs_eta");
59
60 book(m_reco_d0_z0, "reco_d0_vs_z0");
61 book(m_reco_d0_z0sin, "reco_d0_vs_z0sin");
62
63 book(m_truth_pt_vs_eta, "truth_pt_vs_eta");
64 book(m_truth_phi_vs_eta, "truth_phi_vs_eta");
65
66 book(m_truth_hits, "truth_hits");
67 book(m_truth_hits_vs_eta, "truth_hits_vs_eta");
68
69}
70
71void
73
74 // quantities with xAOD::TruthParticle accessors:
75 float eta = particle.eta();
76 float pt = particle.pt() / Gaudi::Units::GeV;
77 float phi = particle.phi();
78
79 static const SG::ConstAccessor<float> d0Acc("d0");
80 static const SG::ConstAccessor<float> z0Acc("z0");
81 static const SG::ConstAccessor<float> thetaAcc("theta");
82 static const SG::ConstAccessor<float> qOverPAcc("qOverP");
83 static const SG::ConstAccessor<float> prodRAcc("prodR");
84 static const SG::ConstAccessor<float> prodZAcc("prodZ");
85 float d0 = d0Acc.isAvailable(particle) ? d0Acc(particle) : -9999.;
86 float z0 = z0Acc.isAvailable(particle) ? z0Acc(particle) : -9999.;
87 float theta = thetaAcc.isAvailable(particle) ? thetaAcc(particle) : -9999.;
88 float z0sin = (z0Acc.isAvailable(particle) && z0Acc(particle)) ? z0 * std::sin(theta) : -9999.;
89 float qOverP = qOverPAcc.isAvailable(particle) ? qOverPAcc(particle) : -9999.;
90 float prodR = prodRAcc.isAvailable(particle) ? prodRAcc(particle) : -9999.;
91 float prodZ = prodZAcc.isAvailable(particle) ? prodZAcc(particle) : -9999.;
92
93 if(d0 > -9000.) fillHisto(m_truth_d0, d0, weight);
94 if(z0 > -9000.) fillHisto(m_truth_z0, z0, weight);
95 if(theta > -9000.) fillHisto(m_truth_theta, theta, weight);
96 if(phi > -9000.) fillHisto(m_truth_phi, phi, weight);
97 if(z0sin > -9000.) fillHisto(m_truth_z0sin, z0sin, weight);
98 if(qOverP > -9000.) fillHisto(m_truth_qoverp, qOverP, weight);
99 if(prodR > -9000.) fillHisto(m_truth_prodR, prodR, weight);
100 if(prodZ > -9000.) fillHisto(m_truth_prodZ, prodZ, weight);
101
102 fillHisto(m_truth_eta, eta, weight);
103 fillHisto(m_truth_pt, pt, weight);
104 fillHisto(m_truth_lowpt, pt, weight);
105
106 fillHisto(m_truth_pt_vs_eta, pt, eta, weight);
108
109 static const SG::AuxElement::ConstAccessor< float > nSilHitsAcc("nSilHits");
110 if (nSilHitsAcc.isAvailable(particle)) {
111 float hits = nSilHitsAcc(particle);
112 fillHisto(m_truth_hits, hits, weight);
113 fillHisto(m_truth_hits_vs_eta, eta, hits, weight);
114 }
115}
116
117void
119
120 float pt = particle.pt() / Gaudi::Units::GeV;
121 float eta = particle.eta();
122 float phi = particle.phi0();
123
124 float chi2 = particle.chiSquared();
125 float ndof = particle.numberDoF();
126 float chi2Overndof = ndof > 0 ? chi2 / ndof : 0;
127
128
129 fillHisto(m_reco_d0, particle.d0(), weight);
130 fillHisto(m_reco_z0, particle.z0(), weight);
131 fillHisto(m_reco_z0sin, particle.z0()* std::sin(particle.theta()), weight);
132
133 fillHisto(m_reco_phi, phi, weight);
134 fillHisto(m_reco_theta, particle.theta(), weight);
135 fillHisto(m_reco_eta, eta, weight);
136 fillHisto(m_reco_qoverp, particle.qOverP(), weight);
137 fillHisto(m_reco_pt, pt, weight);
138 fillHisto(m_reco_lowpt, pt, weight);
139
140 fillHisto(m_reco_pt_vs_eta, pt, eta, weight);
142
143 fillHisto(m_reco_d0_z0, particle.d0(), particle.z0(), weight);
144 fillHisto(m_reco_d0_z0sin, particle.d0(), particle.z0()*std::sin(particle.theta()), weight);
145
146 fillHisto(m_reco_chi2, chi2, weight);
147 fillHisto(m_reco_ndof, ndof, weight);
148 fillHisto(m_reco_chi2Overndof, chi2Overndof, weight);
149
150 std::bitset<xAOD::TrackPatternRecoInfo::NumberOfTrackRecoInfo> patternInfo = particle.patternRecoInfo();
151 for(unsigned int i = 0; i < xAOD::TrackPatternRecoInfo::NumberOfTrackRecoInfo; i++){
152 if(patternInfo.test(i)) fillHisto(m_reco_author, i, weight);
153 }
154
155 if(m_isITk){
156 static const SG::Accessor< uint8_t > accValidTime("hasValidTime");
157 static const SG::Accessor< float > accTime("time");
158 if( accValidTime.isAvailable(particle) && accTime.isAvailable(particle) ) {
159 if (particle.hasValidTime()) {
160 fillHisto(m_reco_time, particle.time(), weight);
161 }
162 }
163
164 if( accValidTime.isAvailable(particle) ) {
165 fillHisto(m_reco_hasValidTime_eff_vs_eta, eta, particle.hasValidTime(), weight);
166 }
167 }
168
169}
170
171void
173
174 //
175 //Pie
176 //
177
178}
Scalar eta() const
pseudorapidity method
Scalar phi() const
phi method
Scalar theta() const
theta method
Helper class to provide constant type-safe access to aux data.
InDetPerfPlot_TrackParameters(InDetPlotBase *pParent, const std::string &dirName, bool isITk)
void fill(const xAOD::TrackParticle &particle, float weight)
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 type-safe access to aux data.
SG::ConstAccessor< T, ALLOC > ConstAccessor
Definition AuxElement.h:569
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.
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.
@ NumberOfTrackRecoInfo
maximum number of enums