ATLAS Offline Software
Loading...
Searching...
No Matches
InDetTrackSelectorTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
6// forward declares
9#include "VxVertex/Vertex.h"
11#include "TrkTrack/Track.h"
12// normal includes
14#include "CLHEP/Matrix/Vector.h"
15
16namespace InDet
17{
18
19//_______________________________________________________________________________
20InDetTrackSelectorTool::InDetTrackSelectorTool(const std::string & t, const std::string & n, const IInterface * p)
21 : AthAlgTool(t,n,p)
22{
23 declareInterface<ITrackSelectorTool>(this);
24}
25
26//_______________________________________________________________________________
28= default;
29
30//_______________________________________________________________________________
32{
34 if (!m_trackSumTool.empty())
35 {
36 if(m_trackSumTool.retrieve().isFailure()) {
37 msg(MSG::ERROR)<<" Unable to retrieve "<<m_trackSumTool<<endmsg;
38 return StatusCode::FAILURE;
39 }
40 else {
42 }
43 }
44
45 if ( m_extrapolator.retrieve().isFailure() ) {
46 msg(MSG::ERROR) << "Failed to retrieve tool " << m_extrapolator << endmsg;
47 return StatusCode::FAILURE;
48 }
49
50 return StatusCode::SUCCESS;
51}
52
53
54//_______________________________________________________________________________
55bool InDetTrackSelectorTool::decision(const Trk::Track & track, const Trk::Vertex * vertex) const
56{
57 // decision based on the track parameters
58 if(!decision(track.perigeeParameters(), vertex, track.info().particleHypothesis()))
59 return false;
60
61 // number of hits, silicon hits, b-layer
62 // first ask track for summary
63 std::unique_ptr<Trk::TrackSummary> summaryUniquePtr;
64 const Trk::TrackSummary * summary = track.trackSummary();
65 if (summary == nullptr && m_trackSumToolAvailable) {
66 summaryUniquePtr = m_trackSumTool->summary(Gaudi::Hive::currentContext(), track);
67 summary = summaryUniquePtr.get();
68 }
69
70 if (nullptr==summary) {
71 ATH_MSG_DEBUG( "Track preselection: cannot create a track summary. This track will not pass." );
72 return false;
73 }
74
75 int nPixelHits = summary->get(Trk::numberOfPixelHits);
76 int nPixelDead = summary->get(Trk::numberOfPixelDeadSensors);
77 if (nPixelDead<0)
78 nPixelDead=0;
79
80 int nInLayerHits = summary->get(Trk::numberOfInnermostPixelLayerHits);
81
82 if(nPixelHits+nPixelDead<m_numberOfPixelHits || nInLayerHits<m_numberOfInLayerHits )
83 return false;
84
85 // all ok
86 return true;
87}
88
89//_______________________________________________________________________________
91{
92 if(!decision(&(track.definingParameters()), vertex, Trk::pion))
93 return false;
94
95 const Trk::TrackSummary * summary = track.trackSummary();
96 if (summary == nullptr) {
97 ATH_MSG_INFO( "TrackParticleBase does not have a Track Summary. Rejected." );
98 return false;
99 }
100 int nPixelHits = summary->get(Trk::numberOfPixelHits);
101 int nPixelDead = summary->get(Trk::numberOfPixelDeadSensors);
102 if (nPixelDead<0)
103 nPixelDead=0;
104
105 int nInLayerHits = summary->get(Trk::numberOfInnermostPixelLayerHits);
106
107 return nPixelHits+nPixelDead>=m_numberOfPixelHits && nInLayerHits>=m_numberOfInLayerHits;
108}
109
110//_______________________________________________________________________________
112{
113 // checking pointer first
114 if(nullptr==track || !track->covariance()) {
115 ATH_MSG_WARNING( "Track preselection: Zero pointer to parameterbase* received (most likely a track without perigee). This track will not pass." );
116 return false;
117 }
118
119 // getting the perigee parameters of the track
120 const Trk::Perigee * perigee(nullptr);
121 if(vertex == nullptr)
122 perigee = dynamic_cast<const Trk::Perigee *>(track);
123 else {
124 Trk::PerigeeSurface perigeeSurface(vertex->position());
125 std::unique_ptr<const Trk::TrackParameters> tmp =
126 m_extrapolator->extrapolate(Gaudi::Hive::currentContext(),
127 *track,
128 perigeeSurface,
130 true,
131 hyp);
132 //release only of right type
133 if (tmp && tmp->associatedSurface().type() == Trk::SurfaceType::Perigee) {
134 perigee = static_cast<const Trk::Perigee*>(tmp.release());
135 }
136 }
137
138 if(nullptr == perigee || !perigee->covariance() ) {
139 ATH_MSG_INFO( "Track preselection: cannot make a measured perigee. This track will not pass." );
140 return false;
141 }
142
143 AmgVector(5) trackParameters = perigee->parameters();
144
145 // d0 and z0 cuts
146 double d0 = trackParameters[Trk::d0];
147 if(std::abs(d0) > m_maxD0) { if(vertex != nullptr) { delete perigee; } return false; }
148
149 double z0 = trackParameters[Trk::z0];
150 if (std::abs(z0)*sin(trackParameters[Trk::theta]) > m_IPz0Max)
151 { if(vertex != nullptr) { delete perigee; } return false; }
152 if (std::abs(z0) > m_maxZ0)
153 { if(vertex != nullptr) { delete perigee; } return false; }
154
155 // transverse momentum
156 double pt = perigee->momentum().perp();
157 if(pt<m_minPt) { if(vertex != nullptr) { delete perigee; } return false; }
158
159 // d0 significance
160 double d0Significance=std::abs(trackParameters[Trk::d0]/sqrt( (*perigee->covariance())(Trk::d0,Trk::d0) ));
161 if (d0Significance>m_maxD0overSigmaD0)
162 { if(vertex != nullptr) { delete perigee; } return false; }
163
164 if(vertex != nullptr) { delete perigee; }
165 return true;
166} //end of selection method
167
168} //end of namespace definitions
#define endmsg
#define ATH_MSG_INFO(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
#define AmgVector(rows)
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
MsgStream & msg() const
virtual bool decision(const Trk::Track &track, const Trk::Vertex *vertex) const override
ToolHandle< Trk::ITrackSummaryTool > m_trackSumTool
InDetTrackSelectorTool(const std::string &t, const std::string &n, const IInterface *p)
ToolHandle< Trk::IExtrapolator > m_extrapolator
virtual StatusCode initialize() override
const Amg::Vector3D & momentum() const
Access method for the momentum.
Class describing the Line to which the Perigee refers to.
A summary of the information contained by a track.
This class is a simplest representation of a vertex candidate.
Primary Vertex Finder.
@ anyDirection
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee
@ theta
Definition ParamDefs.h:66
@ d0
Definition ParamDefs.h:63
@ z0
Definition ParamDefs.h:64
ParticleHypothesis
Enumeration for Particle hypothesis respecting the interaction with material.
ParametersBase< TrackParametersDim, Charged > TrackParameters
@ numberOfPixelHits
number of pixel layers on track with absence of hits
@ numberOfInnermostPixelLayerHits
these are the hits in the 1st pixel layer
@ numberOfPixelDeadSensors
number of pixel hits with broad errors (width/sqrt(12))