ATLAS Offline Software
Loading...
Searching...
No Matches
InDetAdaptiveMultiSecVtxFinderTool.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2019-2023 CERN for the benefit of the ATLAS collaboration
3*/
4// Author: Neza Ribaric <neza.ribaric@cern.ch>
5
6/***************************************************************************
7 InDetAdaptiveMultiSecVtxFinderTool.cxx - Description
8 -------------------
9 begin : 01-12-2022
10 authors : Neza Ribaric ( Lancaster University )
11 information : Tool for Secondary Vertex Finding using AdaptiveMultivertexFitter and InDetTrackSelection
12 ***************************************************************************
13
14 *
15 * This class provides an implementation for a secondary
16 * vertex finding tool, which uses the Adaptive Multi Vertex
17 * Fitter to find and fit multiple secondary vertices.
18 *
19 ***************************************************************************/
20
23#include "GaudiKernel/ServiceHandle.h"
24#include "GaudiKernel/ToolHandle.h"
29#include "TrkParticleBase/TrackParticleBaseCollection.h" // type def ...
30#include "TrkTrack/TrackCollection.h" // type def ...
38
39namespace Trk {
40
41 class Track;
43 class ITrackLink;
47} // namespace Trk
48//
49//
50namespace InDet {
52
54 : public extends<AthAlgTool, IAdaptiveMultiSecVertexFinder> // mutable variables are used without protection.
55 {
56 public:
60 using extends::extends;
61
62 InDetAdaptiveMultiSecVtxFinderTool(const std::string& t, const std::string& n, const IInterface* p);
63
67
69
70 StatusCode initialize() override;
71 StatusCode finalize() override;
72
73 std::pair<xAOD::VertexContainer*, xAOD::VertexAuxContainer*> findVertex(
74 const xAOD::TrackParticleContainer* trackParticles) override;
75
76 void setPrimaryVertexPosition(double, double, double) override;
77
78 private:
79 std::pair<xAOD::VertexContainer*, xAOD::VertexAuxContainer*> doVertexing(const std::vector<Trk::ITrackLink*>& trackVector);
80 float findCompatibleTracks(Amg::Vector3D& seedVertex, Trk::ITrackLink* trkIter) const;
81 void countTracksAndNdf(xAOD::Vertex* myxAODVertex, float& ndf, int& ntracks) const;
82 bool checkFit(xAOD::Vertex* actualCandidate) const;
83 int removeTracksFromSeeds(xAOD::Vertex* actualCandidate, std::vector<Trk::ITrackLink*>& seedTracks) const;
84 void removeClosestTrack(Amg::Vector3D& seedVertex, std::vector<Trk::ITrackLink*>& seedTracks, int& nFound) const;
85
86 bool V0check(const std::vector<Amg::Vector3D>& momenta, const Amg::Vector3D& posi) const;
87 const std::vector<Amg::Vector3D> getVertexMomenta(xAOD::Vertex* myxAODVertex) const;
88
89
90 ToolHandle<Trk::AdaptiveMultiVertexFitter> m_VertexFitter{this, "VertexFitterTool", "Trk::AdaptiveMultiVertexFitter", "Multi Vertex Fitter"};
91 ToolHandle<InDet::IInDetTrackSelectionTool> m_trkFilter{this, "BaseTrackSelector", "InDet::DetailedTrackSelectToolRelax", "Base track selection tool"};
92 ToolHandle<InDet::IInDetTrackSelectionTool> m_SVtrkFilter{this, "SecVtxTrackSelector", "InDet::SecVtxTrackSelector", "SV track selection tool"};
93
94 ToolHandle<Trk::IVertexSeedFinder> m_SeedFinder{this, "SeedFinder", "Trk::IndexedCrossDistancesSeedFinder", "Seed finder"};
95 ToolHandle<Trk::IImpactPoint3dEstimator> m_ImpactPoint3dEstimator{this, "ImpactPoint3dEstimator", "Trk::ImpactPoint3dEstimator", "Impact point estimator"};
96
97 // declareInterface<IAdaptiveMultiSecVertexFinder>(this);
98 FloatProperty m_privtxRef{this, "MomentumProjectionOnDirection", -999999.9, "pri vtx ref"};
99 DoubleProperty m_significanceCutSeeding{this, "significanceCutSeeding", 10, "significanceCutSeeding"};
100 DoubleProperty m_minWghtAtVtx{this, "minTrackWeightAtVtx", 0., "minTrackWeightAtVtx"};
101 DoubleProperty m_maxIterations{this, "maxVertices", 25, "max iterations"};
103
104 }; // end of class definitions
105} // namespace InDet
Define macros for attributes used to control the static checker.
#define ATLAS_NOT_THREAD_SAFE
getNoisyStrip() Find noisy strips from hitmaps and write out into xml/db formats
Interface for track selection tool.
ToolHandle< InDet::IInDetTrackSelectionTool > m_trkFilter
std::pair< xAOD::VertexContainer *, xAOD::VertexAuxContainer * > doVertexing(const std::vector< Trk::ITrackLink * > &trackVector)
float findCompatibleTracks(Amg::Vector3D &seedVertex, Trk::ITrackLink *trkIter) const
ToolHandle< InDet::IInDetTrackSelectionTool > m_SVtrkFilter
ToolHandle< Trk::IImpactPoint3dEstimator > m_ImpactPoint3dEstimator
ToolHandle< Trk::AdaptiveMultiVertexFitter > m_VertexFitter
void countTracksAndNdf(xAOD::Vertex *myxAODVertex, float &ndf, int &ntracks) const
const std::vector< Amg::Vector3D > getVertexMomenta(xAOD::Vertex *myxAODVertex) const
virtual ~InDetAdaptiveMultiSecVtxFinderTool()=default
Destructor.
InDetAdaptiveMultiSecVtxFinderTool(const std::string &t, const std::string &n, const IInterface *p)
int removeTracksFromSeeds(xAOD::Vertex *actualCandidate, std::vector< Trk::ITrackLink * > &seedTracks) const
bool V0check(const std::vector< Amg::Vector3D > &momenta, const Amg::Vector3D &posi) const
void setPrimaryVertexPosition(double, double, double) override
void removeClosestTrack(Amg::Vector3D &seedVertex, std::vector< Trk::ITrackLink * > &seedTracks, int &nFound) const
std::pair< xAOD::VertexContainer *, xAOD::VertexAuxContainer * > findVertex(const xAOD::TrackParticleContainer *trackParticles) override
An Abstract Base Class for the LinearizedTrackFactories.
An abstract base class for implementation of Linearization point finders.
Eigen::Matrix< double, 3, 1 > Vector3D
Primary Vertex Finder.
Ensure that the ATLAS eigen extensions are properly loaded.
Vertex_v1 Vertex
Define the latest version of the vertex class.
TrackParticleContainer_v1 TrackParticleContainer
Definition of the current "TrackParticle container version".
void initialize()