ATLAS Offline Software
Loading...
Searching...
No Matches
JetFitterTrackSelectorTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
5#include <cassert>
6
7using namespace InDet;
8
9 JetFitterTrackSelectorTool::JetFitterTrackSelectorTool(const std::string &t, const std::string &n, const IInterface *p)
10 : AthAlgTool(t, n, p)
11 {
12 declareInterface< JetFitterTrackSelectorTool >(this);
13 }
14
16
18
19 if ( m_trkFilter.retrieve().isFailure() ) {
20 msg(MSG::ERROR) << " Unable to retrieve InDet::InDetDetailedTrackSelectorTool" << endmsg;
21 return StatusCode::FAILURE;
22 }
23
24 if ( m_jetFitterUtils.retrieve().isFailure() ) {
25 msg(MSG::ERROR) << " Unable to retrieve InDet::InDetJetFitterUtils/InDetJetFitterUtils" << endmsg;
26 return StatusCode::FAILURE;
27 }
28
29 if ( m_extrapolator.retrieve().isFailure() ) {
30 msg(MSG::ERROR) << " Unable to retrieve Trk::Extrapolator/InDetExtrapolator" << endmsg;
31 return StatusCode::FAILURE;
32 }
33
34 return StatusCode::SUCCESS;
35 }
36
37
39 const TLorentzVector &jetMomentum,
40 const std::vector<const xAOD::IParticle *> &inputTracks) const {
41 // perform the track selection
42 // step 1, apply a track filter like "InDet::InDetDetailedTrackSelectorTool"
43 // step 2, calculate the compatibility of filtered tracks with primary vertex
44 // use this to deduce primaryTracks and secondaryTracks
45
46 ATH_MSG_DEBUG( "Doing track selection on " << inputTracks.size() << " tracks ... " );
47
48 // We need to use normal pointers instead of smart pointers since the code breaks.
49 // We have to fix this issue in the future
50 // if ( m_selectedTracks != nullptr ) delete m_selectedTracks; // May this break the code?
52
53 // Vectors of Trk::ITrackLink to be given to m_selectedTracks once we understand if they are primary of secondary tracks
54 std::vector< const Trk::ITrackLink* > primaryTrackLinks;
55 std::vector< const Trk::ITrackLink* > secondaryTrackLinks;
56
57 // Running on input tracks
58 std::vector<const xAOD::IParticle *>::const_iterator trk_iter = inputTracks.begin();
59 std::vector<const xAOD::IParticle*>::const_iterator trk_end = inputTracks.end();
60
61 int counter = 0;
62 for ( ; trk_iter != trk_end; ++trk_iter ) {
63 // Convert xAOD::IParticle to xAOD::TrackParticle
64 const xAOD::TrackParticle * tmp = dynamic_cast< const xAOD::TrackParticle* > ( *trk_iter );
65 assert( tmp != nullptr ); // in principle should really check that inputTracks only contains TrackParticle objects
66
67 // Compute compatibility and understand track type
68 // -1: track filter failed
69 // 0: extrapolation of MeasuredPerigee failed
70 // 1: primary
71 // 2: secondary
72 int type = computeTrackCompatibility( primaryVertex,jetMomentum,*tmp );
73
74 // Create Trk::ITrackLink collections to be given to selected tracks
75 if (type==1 || type==2) {
77 linkTP.setElement( tmp );
79
80 if ( type == 1) primaryTrackLinks.push_back( link );
81 else if ( type == 2 ) secondaryTrackLinks.push_back( link );
82 }
83 else {
84 continue;
85 }
86
87 // How many tracks we are selecting
88 counter++;
89 }
90
91 ATH_MSG_DEBUG( " Total of selected tracks: "<< counter );
92
93 selectedTracks->setPrimaryTrackLinks( primaryTrackLinks );
94 selectedTracks->setSecondaryTrackLinks( secondaryTrackLinks );
95 return selectedTracks;
96 }
97
98
101 {
102 std::vector<std::string> out;
103 std::string toolname = this->name();
104 std::string delimiter = "_";
105 std::string::size_type firstDelimiter = toolname.find(delimiter);
106 std::string sub = toolname.substr(0, firstDelimiter);
108 sub += "FLIP_SIGN";
109 }
110 return std::string("JetFitter_TrackCompatibility_") + sub;
111 }
112
116 std::vector<std::string>
118 {
119 return std::vector<std::string> { decorationName() };
120 }
121
123 const TLorentzVector &jetMomentum,
124 const xAOD::TrackParticle &track ) const {
125
126 // Decorators for tracks
127 SG::AuxElement::Decorator< float > compatibilityDecorator(decorationName());
128
129 // Apply track filter
130 if ( !m_trkFilter->decision( track, &primaryVertex ) ) {
131 compatibilityDecorator ( track ) = 0.;
132 return -1;
133 }
134
135 // Recomputing Perigee w.r.t PV
136 Trk::PerigeeSurface mySurface( primaryVertex.position() );
137 std::unique_ptr<const Trk::TrackParameters> myMeasuredPerigee(m_extrapolator->extrapolate(
138 Gaudi::Hive::currentContext(),track.perigeeParameters(),mySurface ));
139
140 if ( !myMeasuredPerigee) {
141 ATH_MSG_DEBUG( " Extrapolation to primary vertex failed. Skipping track " );
142 compatibilityDecorator ( track ) = 0.;
143 return 0;
144 }
145
146
147 // Prepare for using jetFitterUtils (for the computation of the compatibility)
148 // Is this conrvertion really necessary?
149 Trk::RecVertex primaryVertexRecVertex( primaryVertex.position(),
150 primaryVertex.covariancePosition(),
151 primaryVertex.numberDoF(),
152 primaryVertex.chiSquared());
153
154 Amg::Vector3D jetMomSpatial( jetMomentum.X(),jetMomentum.Y(),jetMomentum.Z() );
155 double compatibilityValue = m_jetFitterUtils->compatibility( *myMeasuredPerigee,primaryVertexRecVertex ).first;
156 compatibilityValue = fabs( compatibilityValue ) * InDet::InDetJetFitterUtils::get3DLifetimeSignOfTrack( *myMeasuredPerigee,
157 jetMomSpatial,
158 primaryVertexRecVertex );
159
160 // Decorate
161 ATH_MSG_DEBUG( "compatibilityValue = " << compatibilityValue );
162 compatibilityDecorator ( track ) = compatibilityValue;
163
164 // Understand if primary or secondary track particle
165 double cutCompatibilityPVforPosTracks = m_cutCompatibilityPrimaryVertexForPositiveLifetimeTracks;
166 double cutCompatibilityPVforNegTracks = m_cutCompatibilityPrimaryVertexForNegativeLifetimeTracks;
167
169 cutCompatibilityPVforNegTracks = m_cutCompatibilityPrimaryVertexForPositiveLifetimeTracks;
170 cutCompatibilityPVforPosTracks = m_cutCompatibilityPrimaryVertexForNegativeLifetimeTracks;
171 }
172
173 if ( ( compatibilityValue < 0 &&
174 TMath::Prob( fabs( compatibilityValue ),2 ) < cutCompatibilityPVforNegTracks) ||
175 ( compatibilityValue >= 0 &&
176 TMath::Prob( fabs( compatibilityValue ),2 ) < cutCompatibilityPVforPosTracks ) )
177 return 2;
178 else return 1;
179
180 }
181
182
183
#define endmsg
#define ATH_MSG_DEBUG(x)
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
MsgStream & msg() const
static double get3DLifetimeSignOfTrack(const Trk::TrackParameters &track, const Amg::Vector3D &jetMomentum, const Trk::RecVertex &primaryVertex)
Gaudi::Property< double > m_cutCompatibilityPrimaryVertexForPositiveLifetimeTracks
std::vector< std::string > trackDecorationNames() const
Return a list of the names of track decorations created by this tool, in order to allow them to be lo...
JetFitterTrackSelectorTool(const std::string &t, const std::string &n, const IInterface *p)
std::string decorationName() const
Return the name of the decoration we produce.
ToolHandle< Trk::ITrackSelectorTool > m_trkFilter
Gaudi::Property< bool > m_revertFromPositiveToNegativeTags
Gaudi::Property< double > m_cutCompatibilityPrimaryVertexForNegativeLifetimeTracks
int computeTrackCompatibility(const xAOD::Vertex &primaryVertex, const TLorentzVector &jetMomentum, const xAOD::TrackParticle &trk_iter) const
const Trk::SelectedTracksInJet * doTrackSelection(const xAOD::Vertex &primaryVertex, const TLorentzVector &jetMomentum, const std::vector< const xAOD::IParticle * > &inputTracks) const
ToolHandle< InDet::InDetJetFitterUtils > m_jetFitterUtils
ToolHandle< Trk::IExtrapolator > m_extrapolator
SG::Decorator< T, ALLOC > Decorator
Definition AuxElement.h:575
Element link to XAOD TrackParticle.
Class describing the Line to which the Perigee refers to.
Trk::RecVertex inherits from Trk::Vertex.
Definition RecVertex.h:44
void setPrimaryTrackLinks(std::vector< const ITrackLink * > &primaryTrackLinks)
Set the primary tracks (takes ownership of pointers)
void setSecondaryTrackLinks(std::vector< const ITrackLink * > &secondaryTracLinks)
Set the secondary tracks (takes ownership of pointers)
float numberDoF() const
Returns the number of degrees of freedom of the vertex fit as float.
float chiSquared() const
Returns the of the vertex fit as float.
const Amg::Vector3D & position() const
Returns the 3-pos.
Eigen::Matrix< double, 3, 1 > Vector3D
Primary Vertex Finder.
TrackParticle_v1 TrackParticle
Reference the current persistent version:
Vertex_v1 Vertex
Define the latest version of the vertex class.