ATLAS Offline Software
Loading...
Searching...
No Matches
VertexFittingTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5
6// Local
9
10// Athena
12
13// C/C++
14#include <cmath>
15#include <iostream>
16#include <sstream>
17
18using namespace std;
19
20//=============================================================================
22 const std::string& t, const std::string &name, const IInterface* p
23): AthAlgTool(t, name, p)
24{
25 declareInterface<Prompt::IVertexFittingTool>(this);
26}
27
28//=============================================================================
30{
31 ATH_CHECK(m_vertexFitter.retrieve());
32
34 ATH_CHECK(m_seedVertexFitter.retrieve());
35 }
36
37 m_distToPriVtx = std::make_unique<SG::AuxElement::Decorator<float> > (m_distToPriVtxName);
38 m_normDistToPriVtx = std::make_unique<SG::AuxElement::Decorator<float> > (m_normDistToPriVtxName);
39 m_distToRefittedPriVtx = std::make_unique<SG::AuxElement::Decorator<float> > (m_distToRefittedPriVtxName);
40 m_normDistToRefittedPriVtx = std::make_unique<SG::AuxElement::Decorator<float> > (m_normDistToRefittedPriVtxName);
41 m_distToRefittedRmLepPriVtx = std::make_unique<SG::AuxElement::Decorator<float> > (m_distToRefittedRmLepPriVtxName);
42 m_normDistToRefittedRmLepPriVtx = std::make_unique<SG::AuxElement::Decorator<float> > (m_normDistToRefittedRmLepPriVtxName);
43
44 m_timer.Reset();
45
46 ATH_MSG_DEBUG("VertexFittingSvc::initialize - doSeedVertexFit = " << m_doSeedVertexFit);
47 ATH_MSG_DEBUG("VertexFittingSvc::initialize - vertexFitter = " << m_vertexFitter.name());
48
49 return StatusCode::SUCCESS;
50}
51
52//=============================================================================
54{
55 //
56 // Delete pointers
57 //
58 m_timer.Stop();
59
60 ATH_MSG_INFO("VertexFittingSvc::finalize - number of total fits: " << m_countNumberOfFits);
61 ATH_MSG_INFO("VertexFittingSvc::finalize - number of failed fits: " << m_countNumberOfFitsFailed);
62 ATH_MSG_INFO("VertexFittingSvc::finalize - number of invalid vtxs: " << m_countNumberOfFitsInvalid);
63 ATH_MSG_INFO("VertexFittingSvc::finalize - fitting timer: " << PrintResetStopWatch(m_timer));
64
65 return StatusCode::SUCCESS;
66}
67
68//=============================================================================
70 const FittingInput &input,
71 const std::vector<const xAOD::TrackParticle* > &tracks,
72 VtxType vtxType
73)
74{
75 //
76 // Make new secondary vertex
77 //
80
81 if(!input.priVtx) {
82 ATH_MSG_WARNING("fitVertexWithPrimarySeed -- missing primary vertex");
83 return nullptr;
84 }
85
86 // Fit a new secondary vertex
87 std::unique_ptr<xAOD::Vertex> secondaryVtx = getSecondaryVertexWithSeed(
88 tracks, input.inDetTracks, input.priVtx->position()
89 );
90
91 if(!secondaryVtx) {
93 ATH_MSG_WARNING("fitVertexWithPrimarySeed -- failed to fit vertex");
94 return nullptr;
95 }
96
97 if(!isValidVertex(secondaryVtx.get())) {
99
100 ATH_MSG_WARNING("fitVertexWithPrimarySeed -- failed to get valid vertex");
101 return nullptr;
102 }
103
105
106 // Decorate the newly created vertex
107 static const SG::AuxElement::Accessor<int> indexAcc("SecondaryVertexIndex");
108 static const SG::AuxElement::Accessor<int> typeAcc("SVType");
109
110 indexAcc(*secondaryVtx) = m_secondaryVertexIndex;
111 typeAcc(*secondaryVtx) = static_cast<int>(vtxType);
112
113 decorateNewSecondaryVertex(input, secondaryVtx.get());
114
115 return secondaryVtx;
116}
117
118//=============================================================================
120 const FittingInput &input,
121 const std::vector<const xAOD::TrackParticle* > &tracks,
122 const Amg::Vector3D& seed,
123 VtxType vtxType
124)
125{
126 //
127 // Make new secondary vertex
128 //
131
132 std::unique_ptr<xAOD::Vertex> secondaryVtx = getSecondaryVertexWithSeed(
133 tracks, input.inDetTracks, seed
134 );
135
136 if(!secondaryVtx) {
138
139 ATH_MSG_WARNING("fitVertexWithSeed -- failed to fit vertex");
140 return nullptr;
141 }
142
143 if(!isValidVertex(secondaryVtx.get())) {
145
146 ATH_MSG_WARNING("fitVertexWithSeed -- failed to get valid vertex");
147 return nullptr;
148 }
149
151
152 // Decorate the newly created vertex
153 static const SG::AuxElement::Accessor<int> indexAcc("SecondaryVertexIndex");
154 static const SG::AuxElement::Accessor<int> typeAcc("SVType");
155
156 indexAcc(*secondaryVtx) = m_secondaryVertexIndex;
157 typeAcc(*secondaryVtx) = static_cast<int>(vtxType);
158
159 decorateNewSecondaryVertex(input, secondaryVtx.get());
160
161 return secondaryVtx;
162}
163
164//=============================================================================
166{
167 //
168 // Check if vertex is valid
169 //
170 bool bad_vtx = false;
171
172 if(!vtx){
173 ATH_MSG_WARNING("VertexFittingSvc::Validate_Vertex -- invalid vtx pointer!!!");
174 return false;
175 }
176
177 if(vtx->covariance().empty()) {
178 bad_vtx = true;
179 ATH_MSG_WARNING("VertexFittingSvc::Validate_Vertex -- empty vtx covariance!!!");
180 }
181
182 float chisquared = -9999;
183
184 if(!getVar(vtx, chisquared, "chiSquared")) {
185 bad_vtx = true;
186 ATH_MSG_WARNING("VertexFittingSvc::Validate_Vertex -- not valid vtx chiSquared!!!");
187 }
188
189 float numberdof = -9999;
190
191 if(!getVar(vtx, numberdof, "numberDoF")) {
192 bad_vtx = true;
193 ATH_MSG_WARNING("VertexFittingSvc::Validate_Vertex -- not valid vtx numberDoF!!!");
194 }
195
196 if(std::isnan(vtx->x()) || std::isnan(vtx->y()) || std::isnan(vtx->z())) {
197 bad_vtx = true;
198 ATH_MSG_WARNING("VertexFittingSvc::Validate_Vertex -- vertex coordinate is nan");
199 }
200
201 if(bad_vtx) {
202 ATH_MSG_WARNING("VertexFittingSvc::Validate_Vertex -- bad vertex!!!");
203 ATH_MSG_INFO(printPromptVertexAsStr(vtx, msg(MSG::WARNING)));
204 }
205
206 return !bad_vtx;
207}
208
209//=============================================================================
210void Prompt::VertexFittingTool::removeDoubleEntries(std::vector<const xAOD::TrackParticle*>& tracks)
211{
212 const unsigned nbefore = tracks.size();
213
214 sort(tracks.begin(), tracks.end());
215
216 typename std::vector<const xAOD::TrackParticle*>::iterator TransfEnd = std::unique(tracks.begin(), tracks.end());
217
218 tracks.erase(TransfEnd, tracks.end());
219
220 if(nbefore != tracks.size()) {
221 ATH_MSG_DEBUG("removeDoubleEntries nbefore != tracks.size()): " << nbefore << " != " << tracks.size());
222 }
223}
224
225//=============================================================================
227 const FittingInput &input,
228 xAOD::Vertex *secVtx
229)
230{
231 //
232 // Decorate secondary vertex with all useful information
233 //
234 if(!secVtx) {
235 ATH_MSG_WARNING("decorateNewSecondaryVertex - invalid pointer");
236 return false;
237 }
238
239 //
240 // Decorate secondary vertex with the distance/norm_distance it from Prmary Vertex, Refitted Primary vertex and Refitted Primary Vertex that removed lepton itself.
241 //
242 float distToPriVtx = -1;
243 float normDistToPriVtx = -1;
244 float distToRefittedPriVtx = -1;
245 float normDistToRefittedPriVtx = -1;
246 float distToRefittedRmLepPriVtx = -1;
247 float normDistToRefittedRmLepPriVtx = -1;
248
249 if(input.priVtx) {
250 distToPriVtx = Prompt::getDistance(input.priVtx->position(), secVtx->position());
251 normDistToPriVtx = Prompt::getNormDist(input.priVtx->position(), secVtx->position(), secVtx->covariance(), msg(MSG::WARNING));
252 }
253
254 if(input.refittedPriVtx) {
255 distToRefittedPriVtx = Prompt::getDistance(input.refittedPriVtx->position(), secVtx->position());
256 normDistToRefittedPriVtx = Prompt::getNormDist(input.refittedPriVtx->position(), secVtx->position(), secVtx->covariance(), msg(MSG::WARNING));
257 }
258
259 if(input.refittedPriVtxWithoutLep) {
260 distToRefittedRmLepPriVtx = Prompt::getDistance(input.refittedPriVtxWithoutLep->position(), secVtx->position());
261 normDistToRefittedRmLepPriVtx = Prompt::getNormDist(input.refittedPriVtxWithoutLep->position(), secVtx->position(), secVtx->covariance(), msg(MSG::WARNING));
262 }
263
264 (*m_distToPriVtx) (*secVtx) = distToPriVtx;
265 (*m_normDistToPriVtx) (*secVtx) = normDistToPriVtx;
266 (*m_distToRefittedPriVtx) (*secVtx) = distToRefittedPriVtx;
267 (*m_normDistToRefittedPriVtx) (*secVtx) = normDistToRefittedPriVtx;
268 (*m_distToRefittedRmLepPriVtx) (*secVtx) = distToRefittedRmLepPriVtx;
269 (*m_normDistToRefittedRmLepPriVtx)(*secVtx) = normDistToRefittedRmLepPriVtx;
270
271 return true;
272}
273
274//=============================================================================
276 const std::vector<const xAOD::TrackParticle*> &tracks,
277 const xAOD::TrackParticleContainer *inDetTracks,
278 const Amg::Vector3D& seed
279)
280{
281 //
282 // Fit one vertex with given tracks
283 //
284 std::vector<const xAOD::TrackParticle*> tracksForFit(tracks);
285
286 ATH_MSG_DEBUG("getSecondaryVertexWithSeed -- before remove " << tracksForFit.size());
287
288 removeDoubleEntries(tracksForFit);
289
290 ATH_MSG_DEBUG("getSecondaryVertexWithSeed -- after remove " << tracksForFit.size());
291
292 if(tracksForFit.size() < 2) {
293 ATH_MSG_WARNING("getSecondaryVertexWithSeed -- cannot fit vertex with one or zero input track: ntrack=" << tracksForFit.size());
294 return 0;
295 }
296
297 //
298 // Run fit
299 //
300 ATH_MSG_DEBUG(name() << "::getSecondaryVertexWithSeed -- N tracks = " << tracksForFit.size());
301
302 for(const xAOD::TrackParticle* track: tracksForFit) {
303 ATH_MSG_DEBUG( name() << "::getSecondaryVertexWithSeed -- track pt, eta = " << track->pt() << "," << track->eta());
304 ATH_MSG_DEBUG( name() << "::getSecondaryVertexWithSeed -- track chi2 = " << track->chiSquared());
305 }
306
307 xAOD::Vertex *newVertex = 0;
308 std::unique_ptr<xAOD::Vertex> seedVertex;
309
311 seedVertex = std::unique_ptr<xAOD::Vertex>(m_seedVertexFitter->fit(tracksForFit, seed));
312
313 if(seedVertex.get() && !isValidVertex(seedVertex.get())) {
314 ATH_MSG_DEBUG("getSecondaryVertexWithSeed -- failed to fit seed vertex");
315
316 seedVertex.reset();
317 }
318 }
319
320 if(seedVertex.get()) {
321 newVertex = m_vertexFitter->fit(tracksForFit, seedVertex->position());
322 }
323 else {
324 newVertex = m_vertexFitter->fit(tracksForFit, seed);
325 }
326
327 if(!newVertex) {
328 ATH_MSG_INFO("getSecondaryVertexWithSeed -- failed to fit vertex and fitter returned null xAOD::Vertex pointer");
329 return 0;
330 }
331
332 //
333 // Save vertex tracks
334 //
335 std::vector<ElementLink< xAOD::TrackParticleContainer> > tpLinks;
336
337 if(inDetTracks) {
338 //
339 // Record links to ID tracks if container pointer is provided
340 //
341 for(const xAOD::TrackParticle *selectedtrack: tracksForFit) {
343
344 tpLink.toContainedElement(*inDetTracks, selectedtrack);
345 tpLinks.push_back(tpLink);
346 }
347 }
348
349 newVertex->setTrackParticleLinks(tpLinks);
350
351 TLorentzVector Momentum;
352
353 for(const xAOD::TrackParticle* track: tracksForFit) {
354 Momentum += static_cast<TLorentzVector>(track->p4());
355 }
356
357 xAOD::SecVtxHelper::setVertexMass(newVertex, Momentum.M()); // "mass"
358
359 // newVertex was returned from the vertex fitter as a raw pointer
360 // It looks like Trk::IVertexFitter::fit function expects
361 // memory management to be done by the caller.
362 // Therefore, we take ownership by casting to a unique_ptr
363
364 std::unique_ptr<xAOD::Vertex> returnPtr(newVertex);
365
366 return returnPtr;
367}
368
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_INFO(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
MsgStream & msg() const
std::unique_ptr< SG::AuxElement::Decorator< float > > m_normDistToRefittedRmLepPriVtx
Gaudi::Property< std::string > m_normDistToRefittedPriVtxName
ToolHandle< Trk::IVertexFitter > m_seedVertexFitter
Gaudi::Property< std::string > m_distToPriVtxName
Gaudi::Property< bool > m_doSeedVertexFit
virtual StatusCode initialize() override
std::unique_ptr< xAOD::Vertex > getSecondaryVertexWithSeed(const std::vector< const xAOD::TrackParticle * > &tracks, const xAOD::TrackParticleContainer *inDetTracks, const Amg::Vector3D &seed)
bool decorateNewSecondaryVertex(const FittingInput &input, xAOD::Vertex *secVtx)
Gaudi::Property< std::string > m_distToRefittedRmLepPriVtxName
virtual std::unique_ptr< xAOD::Vertex > fitVertexWithPrimarySeed(const FittingInput &input, const std::vector< const xAOD::TrackParticle * > &tracks, VtxType vtx) override
virtual StatusCode finalize() override
virtual bool isValidVertex(const xAOD::Vertex *vtx) const override
Gaudi::Property< std::string > m_normDistToRefittedRmLepPriVtxName
ToolHandle< Trk::IVertexFitter > m_vertexFitter
std::unique_ptr< SG::AuxElement::Decorator< float > > m_distToPriVtx
virtual std::unique_ptr< xAOD::Vertex > fitVertexWithSeed(const FittingInput &input, const std::vector< const xAOD::TrackParticle * > &tracks, const Amg::Vector3D &seed, VtxType vtxType) override
void removeDoubleEntries(std::vector< const xAOD::TrackParticle * > &tracks)
std::unique_ptr< SG::AuxElement::Decorator< float > > m_normDistToPriVtx
std::unique_ptr< SG::AuxElement::Decorator< float > > m_normDistToRefittedPriVtx
Gaudi::Property< std::string > m_distToRefittedPriVtxName
std::unique_ptr< SG::AuxElement::Decorator< float > > m_distToRefittedRmLepPriVtx
Gaudi::Property< std::string > m_normDistToPriVtxName
std::unique_ptr< SG::AuxElement::Decorator< float > > m_distToRefittedPriVtx
VertexFittingTool(const std::string &t, const std::string &name, const IInterface *p)
SG::Accessor< T, ALLOC > Accessor
Definition AuxElement.h:572
float z() const
Returns the z position.
void setTrackParticleLinks(const TrackParticleLinks_t &trackParticles)
Set all track particle links at once.
float y() const
Returns the y position.
const Amg::Vector3D & position() const
Returns the 3-pos.
float x() const
Returns the x position.
const std::vector< float > & covariance() const
Returns the covariance matrix as a simple vector of values.
Eigen::Matrix< double, 3, 1 > Vector3D
std::string PrintResetStopWatch(TStopwatch &watch)
std::string printPromptVertexAsStr(const xAOD::Vertex *vtx, MsgStream &msg)
bool getVar(T1 &obj, T2 &value, const std::string &var_name)
Definition PromptUtils.h:62
double getNormDist(const Amg::Vector3D &PrimVtx, const Amg::Vector3D &SecVtx, const std::vector< float > &ErrorMatrix, MsgStream &msg)
double getDistance(const xAOD::Vertex *vtx1, const xAOD::Vertex *vtx2)
STL namespace.
DataModel_detail::iterator< DVL > unique(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of unique for DataVector/List.
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.
void setVertexMass(xAOD::Vertex *, float value)
TrackParticle_v1 TrackParticle
Reference the current persistent version:
Vertex_v1 Vertex
Define the latest version of the vertex class.
TrackParticleContainer_v1 TrackParticleContainer
Definition of the current "TrackParticle container version".
MsgStream & msg
Definition testRead.cxx:32