ATLAS Offline Software
Loading...
Searching...
No Matches
InDetIterativePriVxFinderTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3 */
4
5/***************************************************************************
6 InDetIterativePriVxFinderTool.cxx - Description
7 -------------------
8 begin : 23-02-2010
9 authors : Giacinto Piacquadio (CERN PH-ADE-ID)
10 changes :
11 2016-04-26 David Shope <david.richard.shope@cern.ch>
12 EDM Migration to xAOD - from Trk::VxCandidate to xAOD::Vertex
13 findVertex will now always return an xAOD::VertexContainer,
14 even when using a TrackCollection.
15 as input.
16***************************************************************************/
25#include "TrkTrack/Track.h"
27
29#include "VxVertex/RecVertex.h"
31
35
39#include "xAODTracking/Vertex.h"
42
43#include <utility>
44#include <vector>
45
46namespace InDet {
48 const std::string& t,
49 const std::string& n,
50 const IInterface* p)
51 : AthAlgTool(t, n, p)
52{
53 declareInterface<IVertexFinder>(this);
54}
55
57
58StatusCode
60{
62 ATH_MSG_FATAL(" Split vertices cannot be obtained if beam spot constraint "
63 "is true! Change settings...");
64 return StatusCode::FAILURE;
65 }
66
67 /* Get the right vertex fitting tool from ToolSvc */
68 if (m_iVertexFitter.retrieve().isFailure()) {
69 ATH_MSG_FATAL("Failed to retrieve tool " << m_iVertexFitter);
70 return StatusCode::FAILURE;
71 }
72
73 if (m_SeedFinder.retrieve().isFailure()) {
74 ATH_MSG_FATAL("Failed to retrieve tool " << m_SeedFinder);
75 return StatusCode::FAILURE;
76 }
77
78 if (m_LinearizedTrackFactory.retrieve().isFailure()) {
79 ATH_MSG_FATAL("Failed to retrieve tool " << m_LinearizedTrackFactory);
80 return StatusCode::FAILURE;
81 }
82
83 if (m_ImpactPoint3dEstimator.retrieve().isFailure()) {
84 ATH_MSG_FATAL("Failed to retrieve tool " << m_ImpactPoint3dEstimator);
85 return StatusCode::FAILURE;
86 }
87
88 ATH_CHECK(m_beamSpotKey.initialize());
89
90 if (m_trkFilter.retrieve().isFailure()) {
91 ATH_MSG_ERROR(" Unable to retrieve " << m_trkFilter);
92 return StatusCode::FAILURE;
93 }
94
95 // since some parameters special to an inherited class this method
96 // will be overloaded by the inherited class
98
99 return StatusCode::SUCCESS;
100}
101
102namespace {
103struct xAODVertex_pair
104{
105 double first;
106 xAOD::Vertex* second;
107 xAODVertex_pair(double p1, xAOD::Vertex* p2)
108 : first(p1)
109 , second(p2)
110 {}
111 bool operator<(const xAODVertex_pair& other) const
112 {
113 return first > other.first;
114 }
115};
116} // anonymous namespace
117
118std::pair<xAOD::VertexContainer*, xAOD::VertexAuxContainer*>
120 const TrackCollection* trackTES) const
121{
122
124 " Number of input tracks before track selection: " << trackTES->size());
125
127 const InDet::BeamSpotData* beamSpot = *beamSpotHandle;
128
129 std::vector<Trk::ITrackLink*> selectedTracks;
130
131 using TrackDataVecIter = DataVector<Trk::Track>::const_iterator;
132
133 bool selectionPassed;
134 for (TrackDataVecIter itr = (*trackTES).begin(); itr != (*trackTES).end();
135 ++itr) {
136
137 if (m_useBeamConstraint && beamSpot != nullptr) {
138 Trk::RecVertex beamPosition{ beamSpot->beamVtx() };
139 selectionPassed =
140 static_cast<bool>(m_trkFilter->accept(**itr, &beamPosition));
141 } else {
142 Trk::Vertex null(Amg::Vector3D(0, 0, 0));
143 selectionPassed = static_cast<bool>(m_trkFilter->accept(**itr, &null));
144 }
145 if (selectionPassed) {
147 link.setElement(*itr);
148 Trk::LinkToTrack* linkTT = new Trk::LinkToTrack(link);
149 linkTT->setStorableObject(*trackTES);
150 selectedTracks.push_back(linkTT);
151 }
152 }
153
154 ATH_MSG_DEBUG("Of " << trackTES->size() << " tracks " << selectedTracks.size()
155 << " survived the preselection.");
156
157 std::pair<xAOD::VertexContainer*, xAOD::VertexAuxContainer*>
158 returnContainers = findVertex(ctx, selectedTracks);
159
160 return returnContainers;
161}
162
163std::pair<xAOD::VertexContainer*, xAOD::VertexAuxContainer*>
165 const EventContext& ctx,
166 const xAOD::TrackParticleContainer* trackParticles) const
167{
168 ATH_MSG_DEBUG(" Number of input tracks before track selection: "
169 << trackParticles->size());
170
172 const InDet::BeamSpotData* beamSpot = *beamSpotHandle;
173
174 std::vector<Trk::ITrackLink*> selectedTracks;
175
176 using TrackParticleDataVecIter =
178
179 bool selectionPassed;
180 for (TrackParticleDataVecIter itr = trackParticles->begin();
181 itr != trackParticles->end();
182 ++itr) {
183
184 if (m_useBeamConstraint && beamSpot != nullptr) {
185 xAOD::Vertex beamPosition;
186 beamPosition.makePrivateStore();
187 beamPosition.setPosition(beamSpot->beamVtx().position());
188 beamPosition.setCovariancePosition(
189 beamSpot->beamVtx().covariancePosition());
190 selectionPassed =
191 static_cast<bool>(m_trkFilter->accept(**itr, &beamPosition));
192 } else {
193
194 xAOD::Vertex null;
195 null.makePrivateStore();
196 null.setPosition(Amg::Vector3D(0, 0, 0));
197 AmgSymMatrix(3) vertexError;
198 vertexError.setZero();
199 null.setCovariancePosition(vertexError);
200 selectionPassed = static_cast<bool>(m_trkFilter->accept(**itr, &null));
201 }
202
203 if (selectionPassed) {
205 link.setElement(*itr);
208 linkTT->setStorableObject(*trackParticles);
209 selectedTracks.push_back(linkTT);
210 }
211 }
212
213 ATH_MSG_DEBUG("Of " << trackParticles->size() << " tracks "
214 << selectedTracks.size()
215 << " survived the preselection.");
216
217 std::pair<xAOD::VertexContainer*, xAOD::VertexAuxContainer*>
218 returnContainers = findVertex(ctx, selectedTracks);
219
220 return returnContainers;
221}
222
223std::pair<xAOD::VertexContainer*, xAOD::VertexAuxContainer*>
225 const EventContext& ctx,
226 const std::vector<Trk::ITrackLink*>& trackVector) const
227{
229 const InDet::BeamSpotData* beamSpot = *beamSpotHandle;
230 // two things need to be added
231 // 1) the seeding mechanism
232 // 2) the iterative removal of tracks
233
234 std::vector<Trk::ITrackLink*> origTracks = trackVector;
235 std::vector<Trk::ITrackLink*> seedTracks = trackVector;
236
237 // TODO: xAODVertex_pair never seems to be used?
238 std::vector<xAODVertex_pair> myxAODVertices;
239
240 xAOD::VertexContainer* theVertexContainer = new xAOD::VertexContainer;
241 xAOD::VertexAuxContainer* theVertexAuxContainer =
243 theVertexContainer->setStore(theVertexAuxContainer);
244
245 // bail out early with only Dummy vertex if multiplicity cut is applied and
246 // exceeded
247
248 if (m_doMaxTracksCut && (trackVector.size() > m_maxTracks)) {
249 ATH_MSG_WARNING(trackVector.size()
250 << " tracks - exceeds maximum (" << m_maxTracks
251 << "), skipping vertexing and returning only dummy...");
252 xAOD::Vertex* dummyxAODVertex = new xAOD::Vertex;
253 theVertexContainer->push_back(
254 dummyxAODVertex); // have to add vertex to container here first so it can
255 // use its aux store
256 dummyxAODVertex->setPosition(beamSpot->beamVtx().position());
257 dummyxAODVertex->setCovariancePosition(
258 beamSpot->beamVtx().covariancePosition());
259 dummyxAODVertex->vxTrackAtVertex() = std::vector<Trk::VxTrackAtVertex>();
260 dummyxAODVertex->setVertexType(xAOD::VxType::NoVtx);
261 return std::make_pair(theVertexContainer, theVertexAuxContainer);
262 }
263
264 int iterations = 0;
265 unsigned int seedtracknumber = seedTracks.size();
266
267 // used to store seed info
268 Amg::Vector3D actualVertex;
269
270 // prepare iterators for tracks only necessary for seeding
271 std::vector<Trk::ITrackLink*>::iterator seedBegin;
272 std::vector<Trk::ITrackLink*>::iterator seedEnd;
273
274 do {
275
276 seedBegin = seedTracks.begin();
277 seedEnd = seedTracks.end();
278
279 if (seedtracknumber == 0) {
280 ATH_MSG_DEBUG(" New iteration. No tracks available after track selection "
281 "for seeding. No finding done.");
282 break;
283 }
284
285 ATH_MSG_DEBUG("ITERATION NUMBER " << iterations);
286
287 // now find a new SEED
288
289 std::vector<const Trk::TrackParameters*> perigeeList;
290 for (std::vector<Trk::ITrackLink*>::iterator seedtrkAtVtxIter = seedBegin;
291 seedtrkAtVtxIter != seedEnd;
292 ++seedtrkAtVtxIter) {
293 perigeeList.push_back((*seedtrkAtVtxIter)->parameters());
294 }
295
296 xAOD::Vertex theconstraint;
297 if (m_useBeamConstraint && beamSpot != nullptr) {
298 theconstraint =
299 xAOD::Vertex(); // Default constructor creates a private store
300 theconstraint.setPosition(beamSpot->beamVtx().position());
301 theconstraint.setCovariancePosition(
302 beamSpot->beamVtx().covariancePosition());
303 theconstraint.setFitQuality(
304 beamSpot->beamVtx().fitQuality().chiSquared(),
305 beamSpot->beamVtx().fitQuality().doubleNumberDoF());
306 actualVertex = m_SeedFinder->findSeed(perigeeList, &theconstraint);
307 } else {
308 actualVertex = m_SeedFinder->findSeed(perigeeList);
309 Amg::MatrixX looseConstraintCovariance(3, 3);
310 looseConstraintCovariance.setIdentity();
311 looseConstraintCovariance = looseConstraintCovariance * 1e+8;
312 theconstraint = xAOD::Vertex();
313 theconstraint.setPosition(actualVertex);
314 theconstraint.setCovariancePosition(looseConstraintCovariance);
315 theconstraint.setFitQuality(0., -3.);
316 }
317
318 if (msgLvl(MSG::DEBUG)) {
319 msg(MSG::DEBUG) << " seed at x: " << actualVertex.x()
320 << " at y: " << actualVertex.y()
321 << " at z: " << actualVertex.z() << endmsg;
322 }
323
324 if (actualVertex.z() == 0.) {
325 break;
326 }
327
328 // now question (remove tracks which are too far away??? I think so...)
329 std::vector<const Trk::TrackParameters*> perigeesToFit;
330 std::vector<const Trk::TrackParameters*> perigeesToFitSplitVertex;
331
332 int numberOfTracks(perigeeList.size());
333
334 std::vector<const Trk::TrackParameters*>::const_iterator perigeeListBegin =
335 perigeeList.begin();
336 std::vector<const Trk::TrackParameters*>::const_iterator perigeeListEnd =
337 perigeeList.end();
338
339 int counter = 0;
340
341 for (std::vector<const Trk::TrackParameters*>::const_iterator
342 perigeeListIter = perigeeListBegin;
343 perigeeListIter != perigeeListEnd;
344 ++perigeeListIter) {
345 if (numberOfTracks <= 2) {
346 perigeesToFit.push_back(*perigeeListIter);
347 counter += 1;
348 } else if (numberOfTracks <= 4 && !m_createSplitVertices) {
349 perigeesToFit.push_back(*perigeeListIter);
350 counter += 1;
351 } else if (numberOfTracks <= 4 * m_splitVerticesTrkInvFraction &&
353 // few tracks are left, put them into the fit regardless their position!
354 if (counter % m_splitVerticesTrkInvFraction == 0) {
355 perigeesToFit.push_back(*perigeeListIter);
356 counter += 1;
357 } else {
358 perigeesToFitSplitVertex.push_back(*perigeeListIter);
359 counter += 1;
360 }
361 } else { // check first whether it is not too far away!
362 double distance = 0.;
363 try {
364 std::unique_ptr<Trk::PlaneSurface> mySurface =
365 m_ImpactPoint3dEstimator->Estimate3dIP(
366 *perigeeListIter, &actualVertex, distance);
368 msg(MSG::WARNING) << " ImpactPoint3dEstimator failed to find minimum "
369 "distance between track and vertex seed: "
370 << err.p << endmsg;
371 }
372
373 if (distance < 0) {
374 msg(MSG::WARNING)
375 << " Distance between track and seed vtx is negative: " << distance
376 << endmsg;
377 }
378
379 const Trk::TrackParameters* myPerigee = (*perigeeListIter);
380
381 // very approximate error
382 double error = 0.;
383
384 if (myPerigee && myPerigee->covariance()) {
385 error = std::sqrt((*myPerigee->covariance())(Trk::d0, Trk::d0) +
386 (*myPerigee->covariance())(Trk::z0, Trk::z0));
387 } // end of the security check
388
389 if (error == 0.) {
390 msg(MSG::ERROR) << " Error is zero! " << distance << endmsg;
391 error = 1.;
392 }
393
394 if (distance / error < m_significanceCutSeeding) {
395 if (counter % m_splitVerticesTrkInvFraction == 0 ||
397 perigeesToFit.push_back(*perigeeListIter);
398 counter += 1;
399 } else {
400 perigeesToFitSplitVertex.push_back(*perigeeListIter);
401 counter += 1;
402 }
403 }
404 }
405 }
406
407 if (perigeesToFit.empty()) {
408 ATH_MSG_DEBUG(" No good seed found. Exiting search for vertices...");
409 break;
410 }
411
412 // now you have perigeeToFit and perigeeToFitSplitVertices
413 // AND HERE YOU DO THE FIT
414 // to reassign vertices you look ino what is already in myVxCandidate
415 // you do it only ONCE!
416
417 std::unique_ptr<xAOD::Vertex> myxAODVertex;
418 std::unique_ptr<xAOD::Vertex> myxAODSplitVertex;
419
420 if (m_useBeamConstraint && !perigeesToFit.empty()) {
421 myxAODVertex = m_iVertexFitter->fit(ctx, perigeesToFit, theconstraint);
422 } else if (!m_useBeamConstraint && perigeesToFit.size() > 1) {
423 myxAODVertex = m_iVertexFitter->fit(ctx, perigeesToFit);
424 }
425 if (m_createSplitVertices && perigeesToFitSplitVertex.size() > 1) {
426 myxAODSplitVertex = m_iVertexFitter->fit(ctx, perigeesToFitSplitVertex);
427 }
428
429 double ndf = 0.;
430 int ntracks = 0;
431 countTracksAndNdf(myxAODVertex.get(), ndf, ntracks);
432
433 double ndfSplitVertex = 0.;
434 int ntracksSplitVertex = 0;
435 countTracksAndNdf(myxAODSplitVertex.get(), ndfSplitVertex, ntracksSplitVertex);
436
437 bool goodVertex = myxAODVertex != nullptr &&
438 ((!m_useBeamConstraint && ndf > 0 && ntracks >= 2) ||
439 (m_useBeamConstraint && ndf > 3 && ntracks >= 2));
440
441 if (msgLvl(MSG::DEBUG)) {
442 msg(MSG::DEBUG) << " xAOD::Vertex is pointer: " << myxAODVertex
443 << " ndf = " << ndf << " ntracks (with weight>0.01) "
444 << ntracks << " beam constraint "
445 << (m_useBeamConstraint ? "yes" : "no") << endmsg;
446 }
447
448 if (!goodVertex) {
449 removeAllFrom(perigeesToFit, seedTracks);
450 } else {
452 // now you HAVE a good vertex
453 // but you want to add the tracks which you missed...
454
455 int numberOfAddedTracks = 0;
456
457 // iterate on remaining vertices and cross-check if tracks need to be
458 // attached to new vertex
459 xAOD::VertexContainer::iterator vxBegin = theVertexContainer->begin();
460 xAOD::VertexContainer::iterator vxEnd = theVertexContainer->end();
461
462 for (xAOD::VertexContainer::iterator vxIter = vxBegin; vxIter != vxEnd;
463 ++vxIter) {
464 // now iterate on tracks at vertex
465
466 std::vector<Trk::VxTrackAtVertex>* myVxTracksAtVtx =
467 (&(*vxIter)->vxTrackAtVertex());
468
469 if (!myVxTracksAtVtx)
470 continue;
471
472 std::vector<Trk::VxTrackAtVertex>::iterator tracksBegin =
473 myVxTracksAtVtx->begin();
474 std::vector<Trk::VxTrackAtVertex>::iterator tracksEnd =
475 myVxTracksAtVtx->end();
476
477 for (std::vector<Trk::VxTrackAtVertex>::iterator tracksIter =
478 tracksBegin;
479 tracksIter != tracksEnd;) {
480 // only try with tracks which are not too tightly assigned to
481 // another vertex
482 if ((*tracksIter).weight() > 0.01) {
483 ++tracksIter;
484 continue;
485 }
486
487 const Trk::TrackParameters* trackPerigee =
488 (*tracksIter).initialPerigee();
489
490 if (trackPerigee == nullptr) {
491 ATH_MSG_ERROR("Cast to perigee gives null pointer");
492 return std::make_pair(theVertexContainer, theVertexAuxContainer);
493 }
494
495 double chi2_newvtx = compatibility(*trackPerigee, *myxAODVertex);
496 double chi2_oldvtx = compatibility(*trackPerigee, *(*vxIter));
497
498 bool remove = false;
499
500 if (chi2_newvtx < chi2_oldvtx) {
501 ATH_MSG_DEBUG(" Found track of old vertex (chi2= "
502 << chi2_oldvtx
503 << ") more compatible to new one (chi2= "
504 << chi2_newvtx << ")");
505
506 perigeesToFit.push_back(trackPerigee);
507 // but you need to re-add it to the seedTracks too...
508
509 bool isFound = false;
510
511 std::vector<Trk::ITrackLink*>::iterator origBegin =
512 origTracks.begin();
513 std::vector<Trk::ITrackLink*>::iterator origEnd =
514 origTracks.end();
515
516 for (std::vector<Trk::ITrackLink*>::iterator origIter = origBegin;
517 origIter != origEnd;
518 ++origIter) {
519 if ((*origIter)->parameters() == trackPerigee) {
520 isFound = true;
521 seedTracks.push_back(*origIter);
522 break;
523 }
524 }
525 if (!isFound) {
526 ATH_MSG_WARNING(" Cannot find old perigee to re-add back to "
527 "seed tracks... ");
528 }
529
530 numberOfAddedTracks += 1;
531
532 remove = true;
533 }
534
535 if (remove) {
536 // now you have to delete the track from the old vertex...
537 // easy???
538 tracksIter = myVxTracksAtVtx->erase(tracksIter);
539 tracksBegin = myVxTracksAtVtx->begin();
540 tracksEnd = myVxTracksAtVtx->end();
541 } else {
542 ++tracksIter;
543 }
544 } // end of iterating on tracks
545 } // end of iterating on already found vertices in event
546
547 // now you have to delete the previous xAOD::Vertex, do a new fit, then
548 // check if you still have a good vertex
549
550 if (numberOfAddedTracks > 0) {
551 if (m_useBeamConstraint && !perigeesToFit.empty()) {
552 myxAODVertex = m_iVertexFitter->fit(ctx, perigeesToFit, theconstraint);
553 } else if (!m_useBeamConstraint && perigeesToFit.size() > 1) {
554 myxAODVertex = m_iVertexFitter->fit(ctx, perigeesToFit);
555 }
556
557 ndf = 0.;
558 ntracks = 0;
559 countTracksAndNdf(myxAODVertex.get(), ndf, ntracks);
560
561 goodVertex = myxAODVertex != nullptr &&
562 ((!m_useBeamConstraint && ndf > 0 && ntracks >= 2) ||
563 (m_useBeamConstraint && ndf > 3 && ntracks >= 2));
564
565 ATH_MSG_DEBUG(" Refitted xAODVertex is pointer: "
566 << myxAODVertex << " ndf = " << ndf
567 << " ntracks (with weight>0.01) " << ntracks
568 << " beam constraint "
569 << (m_useBeamConstraint ? "yes" : "no"));
570
571 if (!goodVertex) {
572 removeAllFrom(perigeesToFit, seedTracks);
573 ATH_MSG_DEBUG(" Adding tracks resulted in an invalid vertex. "
574 "Should be rare. ");
575 }
576 } // end if tracks were added...
577 } // end reassign tracks from previous vertices and refitting if needed
578
579 // need to re-ask goodVertex since it can be changed in the mean time
580 if (goodVertex) {
581 removeCompatibleTracks(myxAODVertex.get(), perigeesToFit, seedTracks);
582 }
583 } // end else case on if not good Vertex
584
585 bool goodSplitVertex = false;
586
588 goodSplitVertex = myxAODSplitVertex != nullptr && ndfSplitVertex > 0 &&
589 ntracksSplitVertex >= 2;
590
591 ATH_MSG_DEBUG(" xAODSplitVertex is pointer: "
592 << myxAODSplitVertex << " ndf = " << ndfSplitVertex
593 << " ntracks (with weight>0.01) " << ntracksSplitVertex);
594
595 if (!goodSplitVertex) {
596 removeAllFrom(perigeesToFitSplitVertex, seedTracks);
597 } else {
599 myxAODSplitVertex.get(), perigeesToFitSplitVertex, seedTracks);
600
601 } // end else if not good Vertex
602 } // end if create split vertices
603
605 if (goodVertex) {
606 theVertexContainer->push_back(std::move(myxAODVertex));
607 }
608 } else {
609 if (goodVertex) {
610 // type does not seem to be set earlier
611 myxAODVertex->setVertexType(xAOD::VxType::PriVtx);
612 theVertexContainer->push_back(std::move(myxAODVertex));
613 } else {
614 xAOD::Vertex* dummyxAODVertex = new xAOD::Vertex;
615 theVertexContainer->push_back(
616 dummyxAODVertex); // have to add vertex to container here first so it
617 // can use its aux store
618 dummyxAODVertex->setPosition(beamSpot->beamVtx().position());
619 dummyxAODVertex->setCovariancePosition(
620 beamSpot->beamVtx().covariancePosition());
621 dummyxAODVertex->vxTrackAtVertex() =
622 std::vector<Trk::VxTrackAtVertex>();
623 dummyxAODVertex->setVertexType(xAOD::VxType::NoVtx);
624 }
625 if (goodSplitVertex) {
626 // type does not seem to be set earlier
627 myxAODSplitVertex->setVertexType(xAOD::VxType::PriVtx);
628 theVertexContainer->push_back(std::move(myxAODSplitVertex));
629 } else {
630 xAOD::Vertex* dummyxAODVertex = new xAOD::Vertex;
631 theVertexContainer->push_back(
632 dummyxAODVertex); // have to add vertex to container here first so it
633 // can use its aux store
634 dummyxAODVertex->setPosition(beamSpot->beamVtx().position());
635 dummyxAODVertex->setCovariancePosition(
636 beamSpot->beamVtx().covariancePosition());
637 dummyxAODVertex->vxTrackAtVertex() =
638 std::vector<Trk::VxTrackAtVertex>();
639 dummyxAODVertex->setVertexType(xAOD::VxType::NoVtx);
640 }
641 }
642
643 ++iterations;
644 } while (seedTracks.size() > 1 && iterations < m_maxVertices);
645
646 if (iterations == m_maxVertices) {
647 ATH_MSG_DEBUG("Reached maximum iterations, have "<<iterations<<" vertices");
648 }
649 else if (iterations > m_maxVertices) {
650 ATH_MSG_WARNING("Exceeded maximum iterations, have "<<iterations<<" vertices, m_maxVertices = "<<m_maxVertices);
651 }
652
653 //---- add dummy vertex at the end
654 //------------------------------------------------------//
655 //---- if one or more vertices are already there: let dummy have same position
656 // as primary vertex
658 if (!theVertexContainer->empty()) {
659 xAOD::Vertex* primaryVtx = theVertexContainer->front();
660 if (!primaryVtx->vxTrackAtVertex().empty()) {
662 xAOD::Vertex* dummyxAODVertex = new xAOD::Vertex;
663 theVertexContainer->push_back(
664 dummyxAODVertex); // have to add vertex to container here first so it
665 // can use its aux store
666 dummyxAODVertex->setPosition(primaryVtx->position());
667 dummyxAODVertex->setCovariancePosition(
668 primaryVtx->covariancePosition());
669 dummyxAODVertex->vxTrackAtVertex() =
670 std::vector<Trk::VxTrackAtVertex>();
671 dummyxAODVertex->setVertexType(xAOD::VxType::NoVtx);
672 } else {
674 }
675 }
676 //---- if no vertex is there let dummy be at beam spot
677
678 else if (theVertexContainer->empty()) {
679 xAOD::Vertex* dummyxAODVertex = new xAOD::Vertex;
680 theVertexContainer->push_back(
681 dummyxAODVertex); // have to add vertex to container here first so it
682 // can use its aux store
683 dummyxAODVertex->setPosition(beamSpot->beamVtx().position());
684 dummyxAODVertex->setCovariancePosition(
685 beamSpot->beamVtx().covariancePosition());
686 dummyxAODVertex->vxTrackAtVertex() = std::vector<Trk::VxTrackAtVertex>();
687 dummyxAODVertex->setVertexType(xAOD::VxType::NoVtx);
688 }
689 // loop over the pile up to set it as pile up (EXCLUDE first and last
690 // vertex, do not do that in split mode)
691 for (unsigned int i = 0; i < theVertexContainer->size() - 1; i++) {
692
694 " Vtx: " << i << " x= " << (*theVertexContainer)[i]->position().x()
695 << " y= " << (*theVertexContainer)[i]->position().y() << " z= "
696 << (*theVertexContainer)[i]->position().z() << " ntracks= "
697 << (*theVertexContainer)[i]->vxTrackAtVertex().size()
698 << " chi2= " << (*theVertexContainer)[i]->chiSquared()
699 << " ndf = " << (*theVertexContainer)[i]->numberDoF());
700 if (i > 0) {
701 (*theVertexContainer)[i]->setVertexType(xAOD::VxType::PileUp);
702 }
703 }
704 }
705
706 xAOD::VertexContainer::iterator vxBegin = theVertexContainer->begin();
707 xAOD::VertexContainer::iterator vxEnd = theVertexContainer->end();
708
709 // prepare iterators for tracks only necessary for seeding
710 std::vector<Trk::ITrackLink*>::iterator origtrkbegin = origTracks.begin();
711 std::vector<Trk::ITrackLink*>::iterator origtrkend = origTracks.end();
712
713 for (xAOD::VertexContainer::iterator vxIter = vxBegin; vxIter != vxEnd;
714 ++vxIter) {
715 std::vector<Trk::VxTrackAtVertex>* myVxTracksAtVtx =
716 &((*vxIter)->vxTrackAtVertex());
717
718 if (!myVxTracksAtVtx)
719 continue;
720
721 std::vector<Trk::VxTrackAtVertex>::iterator tracksBegin =
722 myVxTracksAtVtx->begin();
723 std::vector<Trk::VxTrackAtVertex>::iterator tracksEnd =
724 myVxTracksAtVtx->end();
725
726 bool found = false;
727
728 for (std::vector<Trk::VxTrackAtVertex>::iterator tracksIter = tracksBegin;
729 tracksIter != tracksEnd;
730 ++tracksIter) {
731
732 // now look for corresponding ITrackLink
733 for (std::vector<Trk::ITrackLink*>::iterator origtrkiter = origtrkbegin;
734 origtrkiter != origtrkend;
735 ++origtrkiter) {
736 if ((*origtrkiter)->parameters() == (*tracksIter).initialPerigee()) {
737 found = true;
738
739 // assigning the input track to the fitted vertex through
740 // VxTrackAtVertices vector
741 (*tracksIter).setOrigTrack(*origtrkiter);
742
743 // See if the trklink is to an xAOD::TrackParticle
744 Trk::LinkToXAODTrackParticle* linkToXAODTP = nullptr;
745 Trk::ITrackLink* tmpLink = (*tracksIter).trackOrParticleLink();
746 if (tmpLink->type() == Trk::ITrackLink::ToxAODTrackParticle) {
747 linkToXAODTP = static_cast<Trk::LinkToXAODTrackParticle*>(tmpLink);
748 }
749
750 // If track is an xAOD::TrackParticle, set directly in xAOD::Vertex
751 if (linkToXAODTP) {
752 (*vxIter)->addTrackAtVertex(*linkToXAODTP, (*tracksIter).weight());
753 }
754
755 origTracks.erase(origtrkiter);
756 origtrkbegin = origTracks.begin();
757 origtrkend = origTracks.end();
758 break;
759 }
760 }
761 if (!found) {
762 ATH_MSG_ERROR(" Cannot find vector element to fix links (step 4)! ");
763 }
764 } // end iterate on tracks at vtx
765 } // end iterate on vertices
766
767 for (std::vector<Trk::ITrackLink*>::iterator origtrkiter = origtrkbegin;
768 origtrkiter != origtrkend;
769 ++origtrkiter) {
770 if ((*origtrkiter) != 0) {
771 delete *origtrkiter;
772 *origtrkiter = 0;
773 }
774 }
775
776 return std::make_pair(theVertexContainer, theVertexAuxContainer);
777}
778
779StatusCode
781{
782 return StatusCode::SUCCESS;
783}
784
785void
787{
788 ATH_MSG_DEBUG("VxPrimary initialize(): Parametersettings "
789 << '\n'
790 << "VertexFitter " << m_iVertexFitter << '\n');
791}
792
793
794double
796 const Trk::TrackParameters& measPerigee,
797 const xAOD::Vertex& vertex) const
798{
799 Trk::LinearizedTrack* myLinearizedTrack =
800 m_LinearizedTrackFactory->linearizedTrack(&measPerigee, vertex.position());
801
802 AmgMatrix(2, 2) weightReduced =
803 myLinearizedTrack->expectedCovarianceAtPCA().block<2, 2>(0, 0);
804
805 AmgMatrix(2, 2) errorVertexReduced =
806 (myLinearizedTrack->positionJacobian() *
807 (vertex.covariancePosition() *
808 myLinearizedTrack->positionJacobian().transpose()))
809 .block<2, 2>(0, 0);
810
811 weightReduced += errorVertexReduced;
812
813 weightReduced = weightReduced.inverse().eval();
814 Amg::Vector2D trackParameters2D =
815 myLinearizedTrack->expectedParametersAtPCA().block<2, 1>(0, 0);
816 double returnValue = trackParameters2D.dot(weightReduced * trackParameters2D);
817
818 delete myLinearizedTrack;
819 myLinearizedTrack = nullptr;
820
821 return returnValue;
822}
823
824void
826 std::vector<const Trk::TrackParameters*>& perigeesToFit,
827 std::vector<Trk::ITrackLink*>& seedTracks) const
828{
829 // remove all perigeesToFit and go on...
830
831 std::vector<Trk::ITrackLink*>::iterator seedBegin = seedTracks.begin();
832 std::vector<Trk::ITrackLink*>::iterator seedEnd = seedTracks.end();
833
834 std::vector<const Trk::TrackParameters*>::iterator perigeesToFitBegin =
835 perigeesToFit.begin();
836 std::vector<const Trk::TrackParameters*>::iterator perigeesToFitEnd =
837 perigeesToFit.end();
838
839 for (std::vector<const Trk::TrackParameters*>::iterator perigeesToFitIter =
840 perigeesToFitBegin;
841 perigeesToFitIter != perigeesToFitEnd;
842 ++perigeesToFitIter) {
843
844 bool found = false;
845
846 for (std::vector<Trk::ITrackLink*>::iterator seedIter = seedTracks.begin();
847 seedIter != seedEnd;
848 ++seedIter) {
849 if ((*seedIter)->parameters() == *perigeesToFitIter) {
850 found = true;
851 seedTracks.erase(seedIter);
852 seedBegin = seedTracks.begin();
853 seedEnd = seedTracks.end();
854 break;
855 }
856 }
857 if (!found) {
859 " Cannot find vector element to delete when removing BAD vertex! ");
860 }
861 } // for cycle
862}
863
864void
866 double& ndf,
867 int& ntracks)
868{
869 if (myxAODVertex) {
870 ndf = myxAODVertex->numberDoF();
871
872 std::vector<Trk::VxTrackAtVertex> myVxTracksAtVtx =
873 myxAODVertex->vxTrackAtVertex();
874
875 std::vector<Trk::VxTrackAtVertex>::iterator tracksBegin =
876 myVxTracksAtVtx.begin();
877 std::vector<Trk::VxTrackAtVertex>::iterator tracksEnd =
878 myVxTracksAtVtx.end();
879
880 for (std::vector<Trk::VxTrackAtVertex>::iterator tracksIter = tracksBegin;
881 tracksIter != tracksEnd;
882 ++tracksIter) {
883 if ((*tracksIter).weight() > 0.01) {
884 ntracks += 1;
885 }
886 }
887 }
888}
889
890void
892 xAOD::Vertex* myxAODVertex,
893 std::vector<const Trk::TrackParameters*>& perigeesToFit,
894 std::vector<Trk::ITrackLink*>& seedTracks) const
895{
896 // now you have your new vertex with its new tracks
897 // now you have to get the compatibility also of all tracks which DIDN'T
898 // BELONG to the vertex!
899 std::vector<Trk::VxTrackAtVertex>* tracksAtVertex =
900 &(myxAODVertex->vxTrackAtVertex());
901
902 std::vector<Trk::VxTrackAtVertex>::const_iterator tracksAtVertexBegin =
903 tracksAtVertex->begin();
904 std::vector<Trk::VxTrackAtVertex>::const_iterator tracksAtVertexEnd =
905 tracksAtVertex->end();
906
907 std::vector<Trk::ITrackLink*>::iterator seedBegin = seedTracks.begin();
908 std::vector<Trk::ITrackLink*>::iterator seedEnd = seedTracks.end();
909
910 std::vector<const Trk::TrackParameters*>::iterator perigeesToFitBegin =
911 perigeesToFit.begin();
912 std::vector<const Trk::TrackParameters*>::iterator perigeesToFitEnd =
913 perigeesToFit.end();
914
915 for (std::vector<Trk::VxTrackAtVertex>::const_iterator tracksAtVertexIter =
916 tracksAtVertexBegin;
917 tracksAtVertexIter != tracksAtVertexEnd;
918 ++tracksAtVertexIter) {
919
920 bool found = false;
921
922 for (std::vector<Trk::ITrackLink*>::iterator seedIter = seedBegin;
923 seedIter != seedEnd;
924 ++seedIter) {
925 if ((*seedIter)->parameters() == (*tracksAtVertexIter).initialPerigee()) {
926 found = true;
927 if ((*tracksAtVertexIter).weight() > 0.01) {
928 seedTracks.erase(seedIter);
929 seedBegin = seedTracks.begin();
930 seedEnd = seedTracks.end();
931 }
932 break;
933 }
934 }
935
936 if (!found) {
937 ATH_MSG_ERROR(" Cannot find vector element to delete (step 1)! ");
938 }
939
940 found = false;
941 for (std::vector<const Trk::TrackParameters*>::iterator perigeesToFitIter =
942 perigeesToFitBegin;
943 perigeesToFitIter != perigeesToFitEnd;
944 ++perigeesToFitIter) {
945 if (*perigeesToFitIter == (*tracksAtVertexIter).initialPerigee()) {
946 found = true;
947 if ((*tracksAtVertexIter).weight() > 0.01) {
948 perigeesToFit.erase(perigeesToFitIter);
949 perigeesToFitBegin = perigeesToFit.begin();
950 perigeesToFitEnd = perigeesToFit.end();
951 }
952 break;
953 }
954 }
955
956 if (!found) {
957 ATH_MSG_ERROR(" Cannot find vector element to delete (step 2)! ");
958 }
959 } // finishing iterating on tracks at vertex
960
961 std::vector<Trk::VxTrackAtVertex>* myVxTracksAtVertex =
962 &(myxAODVertex->vxTrackAtVertex());
963
964 std::vector<Trk::VxTrackAtVertex>::iterator tracksBegin =
965 myVxTracksAtVertex->begin();
966 std::vector<Trk::VxTrackAtVertex>::iterator tracksEnd =
967 myVxTracksAtVertex->end();
968
969 for (std::vector<const Trk::TrackParameters*>::iterator perigeesToFitIter =
970 perigeesToFitBegin;
971 perigeesToFitIter != perigeesToFitEnd;
972 ++perigeesToFitIter) {
973 bool found = false;
974
975 // compute the chi2 with respect to the last fitted vertex!
976 //(easy since track was NOT used in the last vertex fit...)
977
978 const Trk::TrackParameters* myPerigee = (*perigeesToFitIter);
979
980 if (myPerigee == nullptr) {
981 ATH_MSG_ERROR(" nullptr to perigee ");
982 return;
983 }
984
985 double chi2 = compatibility(*myPerigee, *myxAODVertex);
986
987 // check if still sufficiently compatible to previous vertex
988 //(CAN BE VERY LOOSE TO BE CONSERVATIVE AGAINST FAR OUTLIERS)
990 for (std::vector<Trk::ITrackLink*>::iterator seedIter =
991 seedTracks.begin();
992 seedIter != seedEnd;
993 ++seedIter) {
994 if ((*seedIter)->parameters() == *perigeesToFitIter) {
995 found = true;
996 seedTracks.erase(seedIter);
997 seedBegin = seedTracks.begin();
998 seedEnd = seedTracks.end();
999 break;
1000 }
1001 }
1002
1003 if (!found) {
1004 ATH_MSG_ERROR(" Cannot find vector element to delete (step 3)! ");
1005 }
1006 } else {
1007 // look if the track is attached to the vertex. If this is the case you
1008 // should delete the track from the vertex!
1009 for (std::vector<Trk::VxTrackAtVertex>::iterator tracksIter = tracksBegin;
1010 tracksIter != tracksEnd;
1011 ++tracksIter) {
1012 if ((*tracksIter).initialPerigee() == *perigeesToFitIter) {
1013 // delete *tracksIter;
1014 // delete has to happen BEFORE the erase (because the iterator will
1015 // point to the next object in the vector AFTER the erase!)
1016 myVxTracksAtVertex->erase(tracksIter);
1017 tracksBegin = myVxTracksAtVertex->begin();
1018 tracksEnd = myVxTracksAtVertex->end();
1019 found = true;
1020 break;
1021 }
1022 }
1023 }
1024 } // iterate on all perigeesToFit
1025}
1026
1027} // end namespace InDet
#define endmsg
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_FATAL(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
An STL vector of pointers that by default owns its pointed-to elements.
bool operator<(const DataVector< T > &a, const DataVector< T > &b)
Vector ordering relation.
#define AmgSymMatrix(dim)
#define AmgMatrix(rows, cols)
size_t size() const
Number of registered mappings.
DataVector< Trk::Track > TrackCollection
This typedef represents a collection of Trk::Track objects.
#define y
#define x
#define z
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
bool msgLvl(const MSG::Level lvl) const
MsgStream & msg() const
DataModel_detail::const_iterator< DataVector > const_iterator
Standard const_iterator.
Definition DataVector.h:838
value_type push_back(value_type pElem)
Add an element to the end of the collection.
DataModel_detail::iterator< DataVector > iterator
Definition DataVector.h:842
const_iterator end() const noexcept
Return a const_iterator pointing past the end of the collection.
const_iterator begin() const noexcept
Return a const_iterator pointing at the beginning of the collection.
const T * front() const
Access the first element in the collection as an rvalue.
size_type size() const noexcept
Returns the number of elements in the collection.
bool empty() const noexcept
Returns true if the collection is empty.
ToolHandle< InDet::IInDetTrackSelectionTool > m_trkFilter
void printParameterSettings()
Internal method to print the parameters setting.
void removeCompatibleTracks(xAOD::Vertex *myxAODVertex, std::vector< const Trk::TrackParameters * > &perigeesToFit, std::vector< Trk::ITrackLink * > &seedTracks) const
static void countTracksAndNdf(const xAOD::Vertex *myxAODVertex, double &ndf, int &ntracks)
InDetIterativePriVxFinderTool(const std::string &t, const std::string &n, const IInterface *p)
Constructor.
ToolHandle< Trk::IVertexFitter > m_iVertexFitter
virtual ~InDetIterativePriVxFinderTool()
Destructor.
double compatibility(const Trk::TrackParameters &measPerigee, const xAOD::Vertex &vertex) const
void removeAllFrom(std::vector< const Trk::TrackParameters * > &perigeesToFit, std::vector< Trk::ITrackLink * > &seedTracks) const
virtual std::pair< xAOD::VertexContainer *, xAOD::VertexAuxContainer * > findVertex(const EventContext &ctx, const TrackCollection *trackTES) const override
Finding method for Trk::Track collection.
ToolHandle< Trk::IVertexSeedFinder > m_SeedFinder
SG::ReadCondHandleKey< InDet::BeamSpotData > m_beamSpotKey
ToolHandle< Trk::IImpactPoint3dEstimator > m_ImpactPoint3dEstimator
ToolHandle< Trk::IVertexLinearizedTrackFactory > m_LinearizedTrackFactory
AUTO - An Undocumented Tracking Object.
Definition LinkToTrack.h:20
Element link to XAOD TrackParticle.
Trk::RecVertex inherits from Trk::Vertex.
Definition RecVertex.h:44
This class is a simplest representation of a vertex candidate.
void setCovariancePosition(const AmgSymMatrix(3)&covariancePosition)
Sets the vertex covariance matrix.
void setVertexType(VxType::VertexType vType)
Set the type of the vertex.
float numberDoF() const
Returns the number of degrees of freedom of the vertex fit as float.
void setPosition(const Amg::Vector3D &position)
Sets the 3-position.
std::vector< Trk::VxTrackAtVertex > & vxTrackAtVertex()
Non-const access to the VxTrackAtVertex vector.
void setFitQuality(float chiSquared, float numberDoF)
Set the 'Fit Quality' information.
const Amg::Vector3D & position() const
Returns the 3-pos.
double chi2(TH1 *h0, TH1 *h1)
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > MatrixX
Dynamic Matrix - dynamic allocation.
Eigen::Matrix< double, 2, 1 > Vector2D
Eigen::Matrix< double, 3, 1 > Vector3D
Primary Vertex Finder.
@ PriVtx
Primary Vertex.
Definition VertexType.h:27
@ d0
Definition ParamDefs.h:63
@ z0
Definition ParamDefs.h:64
ParametersBase< TrackParametersDim, Charged > TrackParameters
VertexType
Vertex types.
@ PileUp
Pile-up vertex.
@ PriVtx
Primary vertex.
@ NoVtx
Dummy vertex. TrackParticle was not used in vertex fit.
VertexAuxContainer_v1 VertexAuxContainer
Definition of the current jet auxiliary container.
VertexContainer_v1 VertexContainer
Definition of the current "Vertex container version".
Vertex_v1 Vertex
Define the latest version of the vertex class.
TrackParticleContainer_v1 TrackParticleContainer
Definition of the current "TrackParticle container version".