ATLAS Offline Software
Loading...
Searching...
No Matches
JetFitterTwoTrackVtxFinderTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
6
9
14
15using namespace InDet;
16
17JetFitterTwoTrackVtxFinderTool::JetFitterTwoTrackVtxFinderTool(const std::string &t, const std::string &n, const IInterface *p)
18 : AthAlgTool(t, n, p),
19 m_tracksDecorator( "VxTrackAtVertex" )
20{
21 declareInterface< JetFitterTwoTrackVtxFinderTool >(this);
22}
23
25
26
28
29 if ( m_CrossDistancesSeedFinder.retrieve().isFailure() ) {
30 ATH_MSG_ERROR( "Cannot retrieve Trk::CrossDistancesSeedFinder/CrossDistancesSeedFinder" );
31 return StatusCode::FAILURE;
32 }
33
34 if ( m_SequentialVertexFitter.retrieve().isFailure() ) {
35 ATH_MSG_ERROR( "Cannot retrieve Trk::SequentialVertexFitter/SequentialVertexFitter" );
36 return StatusCode::FAILURE;
37 }
38
39 return StatusCode::SUCCESS;
40}
41
43 return StatusCode::SUCCESS;
44}
45
46
48 const xAOD::Vertex& primaryVertex,
49 const TLorentzVector& jetMomentum,
50 std::vector< const Trk::ITrackLink* >& inputTracks) const {
51
52 std::vector< const xAOD::Vertex* > VtxCandidates;
53
54 ATH_MSG_DEBUG( "Looping over " << inputTracks.size() << " input tracks... " );
55
56 // Loop over all the combinations
57 for ( unsigned int indexA(0); indexA<inputTracks.size(); indexA++ ) {
58 const Trk::ITrackLink* trackA = inputTracks.at( indexA );
59
60 for ( unsigned int indexB(0); indexB<indexA; indexB++ ) {
61 const Trk::ITrackLink* trackB = inputTracks.at( indexB );
62
63 // Computing the Vertex candidate
64 xAOD::Vertex *myCandidate = computeVtxcandidate( ctx,primaryVertex,jetMomentum,trackA,trackB );
65 if ( myCandidate == nullptr ) continue;
66
67 // Attaching tracks to vertex candidate
68 std::vector< const Trk::ITrackLink* > associatedTracksAtVertex;
69 associatedTracksAtVertex.push_back( trackA );
70 associatedTracksAtVertex.push_back( trackB );
71 m_tracksDecorator ( *myCandidate ) = associatedTracksAtVertex;
72
73 VtxCandidates.push_back( myCandidate );
74 }
75 }
76
77
78 ATH_MSG_DEBUG( "Found " << VtxCandidates.size() <<" 2-trk vertex candidates!" );
79
80 const Trk::TwoTrackVerticesInJet *twoTrackVerticesInJet = new Trk::TwoTrackVerticesInJet ( VtxCandidates,
81 std::vector< const Trk::TrackParticleBase* >() );
82 return twoTrackVerticesInJet;
83}
84
86 const xAOD::Vertex& primaryVertex,
87 const TLorentzVector& jetMomentum,
88 const Trk::ITrackLink* trackA,
89 const Trk::ITrackLink* trackB ) const {
90
91 const Trk::TrackParameters* perigeeTrackA = trackA->parameters();
92 const Trk::TrackParameters* perigeeTrackB = trackB->parameters();
93
94 std::vector< const Trk::TrackParameters* > perigeeToFit;
95 perigeeToFit.push_back( perigeeTrackA );
96 perigeeToFit.push_back( perigeeTrackB );
97
98 Amg::Vector3D seedVertex;
99 try {
100 seedVertex = m_CrossDistancesSeedFinder->findSeed( perigeeToFit );
101
102 if ( seedVertex.perp() > m_maxR ||
103 fabs( seedVertex.z() ) > m_maxZ ) {
104 ATH_MSG_DEBUG( "Vertex seed outside ID. R=" << seedVertex.perp() << " Z=" << seedVertex.z() );
105 seedVertex = primaryVertex.position();
106 }
107
108 } catch (...) {
109 ATH_MSG_DEBUG( "Seed finding failed. Using primary vertex as seed... (clearly not optimal)" );
110 seedVertex = Amg::Vector3D( primaryVertex.position() );
111 }
112
113 ATH_MSG_DEBUG( "Seed: x=" << seedVertex.x() <<
114 " y=" << seedVertex.y() <<
115 " z=" << seedVertex.z() );
116
117
118 // Compute V0 candidate
119 std::unique_ptr< xAOD::Vertex > myCandidate( m_SequentialVertexFitter->fit(ctx,perigeeToFit,seedVertex) );
120
121 // Check fit completed with success
122 if ( myCandidate == nullptr ) {
123 ATH_MSG_DEBUG( " Sequential fit failed. shouldn't happen... Skipping V0 candidate... " );
124 return nullptr;
125 }
126
127 // Check ChiSquared and nDof
128 if ( myCandidate->chiSquared() < 0 ||
129 myCandidate->numberDoF() < 0 ) {
130 ATH_MSG_DEBUG( " Fit for V0 candidate failed: chi2 or ndf negative. Deleting candidate..." );
131 return nullptr;
132 }
133
134 // Check two-vertex probability
135 if ( TMath::Prob( myCandidate->chiSquared(),myCandidate->numberDoF() ) <= m_twoVertexProbabilityCut ) {
136 ATH_MSG_DEBUG( " Candidate not satisfying two-vertex probability " );
137 return nullptr;
138 }
139
140 // Check sign and revertFromPositiveToNegativeTags
141 Amg::Vector3D jetMomSpatial( jetMomentum.X(),jetMomentum.Y(),jetMomentum.Z() );
142 Amg::Vector3D twoTrkVtxPos( myCandidate->position() );
143 double sign = ( twoTrkVtxPos-primaryVertex.position() ).dot( jetMomSpatial );
144
146 ATH_MSG_DEBUG( "Not satisfying sign and revertFromPositiveToNegativeTags requirements" );
147 return nullptr;
148 }
149
150 ATH_MSG_DEBUG( " Candidate: x=" << myCandidate->x() <<
151 " y=" << myCandidate->y() <<
152 " z=" << myCandidate->z() );
153
154 return myCandidate.release();
155}
156
157
#define ATH_MSG_ERROR(x)
#define ATH_MSG_DEBUG(x)
int sign(int a)
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
ToolHandle< Trk::IVertexSeedFinder > m_CrossDistancesSeedFinder
JetFitterTwoTrackVtxFinderTool(const std::string &t, const std::string &n, const IInterface *p)
ToolHandle< Trk::IVertexFitter > m_SequentialVertexFitter
const Trk::TwoTrackVerticesInJet * doVertexFinding(const EventContext &ctx, const xAOD::Vertex &, const TLorentzVector &, std::vector< const Trk::ITrackLink * > &) const
SG::AuxElement::Decorator< std::vector< const Trk::ITrackLink * > > m_tracksDecorator
xAOD::Vertex * computeVtxcandidate(const EventContext &ctx, const xAOD::Vertex &, const TLorentzVector &, const Trk::ITrackLink *trackA, const Trk::ITrackLink *trackB) const
float z() const
Returns the z position.
float y() const
Returns the y position.
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.
float x() const
Returns the x position.
Eigen::Matrix< double, 3, 1 > Vector3D
Primary Vertex Finder.
ParametersBase< TrackParametersDim, Charged > TrackParameters
Vertex_v1 Vertex
Define the latest version of the vertex class.