ATLAS Offline Software
Loading...
Searching...
No Matches
MuSAVtxFitterTool.h
Go to the documentation of this file.
1/*
2Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5#ifndef MUSAVTXFITTERTOOL_H
6#define MUSAVTXFITTERTOOL_H
7
9#include "GaudiKernel/ToolHandle.h"
19
20namespace Rec {
21 static const InterfaceID IID_MuSAVtxFitterTool("MuSAVtxFitterTool", 1, 0);
22
24 {
25 public:
26 MuSAVtxFitterTool(const std::string& type, const std::string& name, const IInterface* parent);
28 virtual StatusCode initialize() override;
29
30 static const InterfaceID& interfaceID() { return IID_MuSAVtxFitterTool;}
31
46
47 struct WrkVrt {
48
49 bool isGood = true;
50 std::deque<long int> selectedTrackIndices;
51 std::deque<long int> associatedTrackIndices;
52 std::vector<std::shared_ptr<const xAOD::TrackParticle>> newExtrapolatedTracks;
53 std::vector<const xAOD::Muon*> muonCandidates;
55 TLorentzVector mom;
56 std::vector<double> cov;
57 double chi2 = 0;
58 double chi2Core = 0;
59 std::vector<double> chi2PerTrk;
60 long int charge = 0;
61 std::vector< std::vector<double> > trkAtVrt;
62 double minOpAng = 0;
64 double closestWrkVrtValue = 0;
65
66 inline double ndof() const { return 2.0*( selectedTrackIndices.size() + associatedTrackIndices.size() ) - 3.0; } //currently always 1, because only 2trk fit, but left for future developments.
67 inline double ndofCore() const { return 2.0*( selectedTrackIndices.size() ) - 3.0; } //currently always 1, because only 2trk fit, but left for future developments. Uses only seeding MSTPs ("core")
68 inline unsigned nTracksTotal() const { return selectedTrackIndices.size() + associatedTrackIndices.size(); }
69 inline double fitQuality() const { return chi2 / ndof(); }
70 };
71
72 // takes in MSTPs from SA muons and extrapolates them backwards without any constraint, producing a new xAOD::TrackParticle for input to the vertex fitter
73 std::unique_ptr<xAOD::TrackParticle> extrapolateMuSA(const xAOD::TrackParticle& trk, const xAOD::EventInfo& eventInfo, const EventContext& ctx) const;
74
75 // performs the MuSA vertex fit routine using the VKalVrtFitter tool
76 StatusCode doMuSAVtxFit(std::vector<MuSAVtxFitterTool::WrkVrt>& workVerticesContainer, const xAOD::MuonContainer& muonContainer, const xAOD::EventInfo& eventInfo, const EventContext& ctx) const;
77
78 // does post-fit selections on vertices in the workVerticesContainer
79 StatusCode doPostFitSelections(std::vector<MuSAVtxFitterTool::WrkVrt>& workVerticesContainer) const;
80
81 // selects the best remaining vertices from the workVerticesContainer and resolves ambiguities to produce unique track-to-vertex associations
82 StatusCode selectBestVertices(std::vector<MuSAVtxFitterTool::WrkVrt>& workVerticesContainer) const;
83
84 private:
85
86 ToolHandle<Trk::ITrkVKalVrtFitter> m_vertexFitter{
87 this, "VertexFitterTool", "Trk::TrkVKalVrtFitter", "Vertex fitter tool"};
88 ToolHandle <Trk::IExtrapolator> m_extrapolator{
89 this,"Extrapolator","Trk::Extrapolator/AtlasExtrapolator","ATLAS Extrapolator tool"};
90
91 DoubleProperty m_etaCutMSTP{this, "etaCutMSTP", 2.5, "Maximum |eta| of input MSTPs"};
92 DoubleProperty m_baseChi2Cut{this, "baseChi2Cut", 50., "Maximum allowed chi2 of saved vertices"};
93 Gaudi::Property<bool> m_doValidation{this, "doValidation", false, "Vertex every MSTP in input as part of validation process"};
94
95 };
96}
97
98#endif
xAOD::MuonContainer * muonContainer
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
StatusCode selectBestVertices(std::vector< MuSAVtxFitterTool::WrkVrt > &workVerticesContainer) const
ToolHandle< Trk::ITrkVKalVrtFitter > m_vertexFitter
virtual StatusCode initialize() override
StatusCode doPostFitSelections(std::vector< MuSAVtxFitterTool::WrkVrt > &workVerticesContainer) const
StatusCode doMuSAVtxFit(std::vector< MuSAVtxFitterTool::WrkVrt > &workVerticesContainer, const xAOD::MuonContainer &muonContainer, const xAOD::EventInfo &eventInfo, const EventContext &ctx) const
static const InterfaceID & interfaceID()
std::unique_ptr< xAOD::TrackParticle > extrapolateMuSA(const xAOD::TrackParticle &trk, const xAOD::EventInfo &eventInfo, const EventContext &ctx) const
Gaudi::Property< bool > m_doValidation
ToolHandle< Trk::IExtrapolator > m_extrapolator
MuSAVtxFitterTool(const std::string &type, const std::string &name, const IInterface *parent)
Eigen::Matrix< double, 3, 1 > Vector3D
Gaudi Tools.
static const InterfaceID IID_MuSAVtxFitterTool("MuSAVtxFitterTool", 1, 0)
EventInfo_v1 EventInfo
Definition of the latest event info version.
TrackParticle_v1 TrackParticle
Reference the current persistent version:
MuonContainer_v1 MuonContainer
Definition of the current "Muon container version".
Represents a working vertex structure used in VKalVrt fitting.
int closestWrkVrtIndex
minimum opening angle between tracks
double chi2
VKalVrt fit covariance.
std::vector< double > cov
VKalVrt fit vertex 4-momentum.
double closestWrkVrtValue
stores the index of the closest WrkVrt in std::vector<WrkVrt> – CURRENTLY UNUSED, LEFT FOR FUTURE DEV...
std::deque< long int > selectedTrackIndices
flagged true for good vertex candidates
long int charge
list of VKalVrt fit chi2 for each track
double minOpAng
list of track parameters wrt the reconstructed vertex
std::deque< long int > associatedTrackIndices
list if indices in TrackParticleContainer for selectedBaseTracks – CURRENTLY UNUSED,...
Amg::Vector3D pos
list of associated muon candidates
std::vector< std::vector< double > > trkAtVrt
total charge of the vertex
std::vector< double > chi2PerTrk
VKalVrt fit chi2 result.
TLorentzVector mom
VKalVrt fit vertex position.
std::vector< const xAOD::Muon * > muonCandidates
list of associated new extrapolated tracks
std::vector< std::shared_ptr< const xAOD::TrackParticle > > newExtrapolatedTracks
list if indices in TrackParticleContainer for associatedTracks – CURRENTLY UNUSED,...
double chi2Core
VKalVrt fit chi2 result.
double ndof() const
stores the value of some observable to the closest WrkVrt ( observable = e.g. significance ) – CURREN...