ATLAS Offline Software
Loading...
Searching...
No Matches
MuSAVtxFitterTool.cxx
Go to the documentation of this file.
1/*
2Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
8
9// Header include
10
13#include "GeoPrimitives/GeoPrimitives.h" //Needed for Amg::Vector3D
14#include <vector>
15#include <cmath>
16
17namespace Rec {
18
19MuSAVtxFitterTool::MuSAVtxFitterTool(const std::string& type, const std::string& name, const IInterface* parent) :
20 AthAlgTool(type, name, parent)
21{
22 declareInterface<MuSAVtxFitterTool>(this);
23}
24
26
28{
29
30 ATH_CHECK(m_vertexFitter.retrieve());
31 ATH_MSG_DEBUG("Retrieved tool " << m_vertexFitter);
32
33 ATH_CHECK(m_extrapolator.retrieve());
34 ATH_MSG_DEBUG("Retrieved tool " << m_extrapolator);
35
36
37 ATH_MSG_DEBUG("Initialization successful");
38
39 return StatusCode::SUCCESS;
40
41}
42
43std::unique_ptr<xAOD::TrackParticle> MuSAVtxFitterTool::extrapolateMuSA(const xAOD::TrackParticle& trk, const xAOD::EventInfo& eventInfo, const EventContext& ctx) const {
44
45 Amg::Vector3D refPos(eventInfo.beamPosX(), eventInfo.beamPosY(), eventInfo.beamPosZ());
46 Trk::PerigeeSurface persurf(refPos);
47 Trk::PropDirection propagationDirection = Trk::oppositeMomentum;
48
49 auto newTrack = std::make_unique<xAOD::TrackParticle>();
50 newTrack->makePrivateStore();
51
52 auto perigee = trk.perigeeParameters();
53 auto entryPars = m_extrapolator->extrapolate(ctx, perigee, persurf, propagationDirection, false, Trk::muon);
54
55 if (!entryPars) {
56 ATH_MSG_WARNING("Failed to back extrapolate MSTP, returning dummy!");
57 // Return a dummy track that will fail downstream cuts (extrapolations can sometimes fail due to various non-pathological reasons -- SA muons are not always high quality objects)
58 newTrack->setDefiningParameters(1.e19, 1.e19, 1.e19, 1.e19, 0);
59 return newTrack;
60 }
61
62 newTrack->setDefiningParameters(entryPars->parameters()[Trk::d0], entryPars->parameters()[Trk::z0], entryPars->parameters()[Trk::phi], entryPars->parameters()[Trk::theta], entryPars->parameters()[Trk::qOverP]);
63 newTrack->setDefiningParametersCovMatrix(*entryPars->covariance());
64 newTrack->setParametersOrigin(eventInfo.beamPosX(), eventInfo.beamPosY(), eventInfo.beamPosZ());
65
66 return newTrack;
67}
68
69StatusCode MuSAVtxFitterTool::doMuSAVtxFit(std::vector<MuSAVtxFitterTool::WrkVrt>& workVerticesContainer,
71 const xAOD::EventInfo& eventInfo,
72 const EventContext& ctx) const
73{
74 ATH_MSG_DEBUG("MuSAVtxFitterTool::doMuSAVtxFit");
75
76 // first gather all SA muons that pass basic checks
77 std::vector<const xAOD::Muon*> candidateSAmuons;
78 for (const auto muon : muonContainer) {
79 bool isSA = (muon->muonType() == xAOD::Muon::MuonStandAlone);
80 bool isCalo = (muon->muonType() == xAOD::Muon::CaloTagged);
81 bool isSegment = (muon->muonType() == xAOD::Muon::SegmentTagged);
82 bool isSiForward = (muon->muonType() == xAOD::Muon::SiliconAssociatedForwardMuon);
83
84 if (!isSA && !m_doValidation) {
85 continue;
86 }
87
88 const xAOD::TrackParticle* MuSAMSTP = muon->trackParticle(xAOD::Muon::MuonSpectrometerTrackParticle);
89 if (!MuSAMSTP) {
90 if (isCalo || isSegment || isSiForward) {
91 ATH_MSG_VERBOSE("Skipping non-SA, non-Combined muon type in validation mode!");
92 } else {
93 ATH_MSG_WARNING("Muon has no MSTP, check your input!");
94 }
95 continue;
96 }
97
98 // pre-fit cuts at "seeding" level go here
99 if (std::abs(MuSAMSTP->eta()) > m_etaCutMSTP) {
100 ATH_MSG_VERBOSE("Skipping SA muon with |eta| > 2.5!");
101 continue;
102 }
103 // SA muons can get saved with nonphysical pT values, checking prevents extrap crashes
104 if (MuSAMSTP->pt() > 13000000) {
105 ATH_MSG_DEBUG("Skipping SA muon with pT " << (MuSAMSTP->pt() / 1000.) << " GeV!");
106 continue;
107 }
108
109 // SA muons can also be saved in regions with 0 magnetic field, which can cause extrapolation crashes
110 float spectrometerFieldIntegral = 0.0;
111 muon->parameter(spectrometerFieldIntegral, xAOD::Muon::spectrometerFieldIntegral);
112 if (spectrometerFieldIntegral < 0.1) {
113 ATH_MSG_DEBUG("Skipping SA muon with spectrometerFieldIntegral " << muon->spectrometerFieldIntegral << " T*m!");
114 continue;
115 }
116
117 candidateSAmuons.push_back(muon);
118 }
119
120 // bail out if we don't have at least two valid SA muons
121 if (candidateSAmuons.size() < 2) {
122 ATH_MSG_VERBOSE("Not enough standalone muons to fit a vertex!");
123 return StatusCode::SUCCESS;
124 }
125
126 // extrapolate SA muons only if we have enough candidates to form a vertex
127 std::vector<std::unique_ptr<xAOD::TrackParticle>> extrapolatedMuSATracks;
128 std::map<const xAOD::TrackParticle*, const xAOD::Muon*> muonToExtrapolatedTrackMap;
129 for (const auto muon : candidateSAmuons) {
130 const xAOD::TrackParticle* MuSAMSTP = muon->trackParticle(xAOD::Muon::MuonSpectrometerTrackParticle);
131 auto extrapolatedMuSATrack = extrapolateMuSA(*MuSAMSTP, eventInfo, ctx);
132 if (extrapolatedMuSATrack->definingParameters()[Trk::d0] > 8000) {
133 ATH_MSG_DEBUG("Failed to extrapolate MuSA track, skipping!");
134 continue;
135 }
136 extrapolatedMuSATracks.push_back(std::move(extrapolatedMuSATrack));
137 muonToExtrapolatedTrackMap[extrapolatedMuSATracks.back().get()] = muon;
138 ATH_MSG_VERBOSE("Extrapolated MuSA track! Total extrapolated so far: " << extrapolatedMuSATracks.size());
139 }
140
141 if (extrapolatedMuSATracks.size() < 2) {
142 ATH_MSG_VERBOSE("Not enough extrapolated tracks to fit a vertex!");
143 return StatusCode::SUCCESS;
144 }
145
146 std::unique_ptr<Trk::IVKalState> state = m_vertexFitter->makeState(ctx);
147 // Loop over all unique pairs.
148 for (unsigned int i = 0; i < extrapolatedMuSATracks.size(); i++) {
149 for (unsigned int j = i+1; j < extrapolatedMuSATracks.size(); j++) {
150 MuSAVtxFitterTool::WrkVrt MuSACandidate;
151 std::vector<const xAOD::TrackParticle*> tracksToFit = {
152 extrapolatedMuSATracks[i].get(),
153 extrapolatedMuSATracks[j].get()
154 };
155 std::vector<const xAOD::NeutralParticle*> dummyNeutrals;
156 StatusCode res = m_vertexFitter->VKalVrtFit(tracksToFit, dummyNeutrals,
157 MuSACandidate.pos, MuSACandidate.mom,
158 MuSACandidate.charge, MuSACandidate.cov,
159 MuSACandidate.chi2PerTrk, MuSACandidate.trkAtVrt,
160 MuSACandidate.chi2, *state);
161 if (res.isSuccess()) {
162 ATH_MSG_DEBUG("MuSA vertex fit successful!");
163 // Transfer persistent ownership by converting the unique_ptr to a shared_ptr
164 // for both tracks used in this candidate (lets unused tracks go out of scope peacefully while keeping used ones)
165 auto sharedTrack1 = std::make_shared<const xAOD::TrackParticle>(*extrapolatedMuSATracks[i]);
166 auto sharedTrack2 = std::make_shared<const xAOD::TrackParticle>(*extrapolatedMuSATracks[j]);
167 MuSACandidate.newExtrapolatedTracks.push_back(sharedTrack1);
168 MuSACandidate.newExtrapolatedTracks.push_back(sharedTrack2);
169 MuSACandidate.muonCandidates = {
170 muonToExtrapolatedTrackMap[tracksToFit[0]],
171 muonToExtrapolatedTrackMap[tracksToFit[1]]
172 };
173 MuSACandidate.minOpAng = muonToExtrapolatedTrackMap[tracksToFit[0]]->p4().DeltaR(
174 muonToExtrapolatedTrackMap[tracksToFit[1]]->p4());
175 workVerticesContainer.emplace_back(MuSACandidate);
176 } else {
177 ATH_MSG_DEBUG("MuSA vertex fit failed!");
178 }
179 }
180 }
181 // NOTE: Any remaining extrapolated tracks in the vector will be cleaned up automatically when they go out of scope
182 // this means the logic may need to be revisited if we want to "recover" extrapolated tracks that were not used in any fit through another mechanism
183 return StatusCode::SUCCESS;
184}
185
186StatusCode MuSAVtxFitterTool::doPostFitSelections(std::vector<MuSAVtxFitterTool::WrkVrt>& workVerticesContainer) const{
187 ATH_MSG_DEBUG("MuSAVtxFitterTool::doPostFitSelections");
188 //minimal method for post-fit selections
189 for (auto& workVertex : workVerticesContainer) {
190 if (workVertex.chi2 > m_baseChi2Cut) {
191 ATH_MSG_DEBUG("SA vertex chi2 too high, removing from consideration");
192 workVertex.isGood = false; // isGood flag has no specific utility here yet, but may be used in the future for more sophisticated recovery/reassembly
193 }
194 }
195
196 //remove any vertices marked as bad
197 workVerticesContainer.erase(std::remove_if(workVerticesContainer.begin(), workVerticesContainer.end(),
198 [&](const MuSAVtxFitterTool::WrkVrt& vtx) { return !vtx.isGood; }), workVerticesContainer.end());
199
200 return StatusCode::SUCCESS;
201}
202
203StatusCode MuSAVtxFitterTool::selectBestVertices(std::vector<MuSAVtxFitterTool::WrkVrt>& workVerticesContainer) const
204{
205 ATH_MSG_DEBUG("MuSAVtxFitterTool::selectBestVertices");
206 //the default logic will create vertices with all possible combinations of SA muons
207 //in rare cases this will cause two MuSA vertices in one event that re-use a muon (more likely in signal models with two or more LLP BSM decays)
208 //if this occurs we can select based on the highest chi2 and remove the poorer-quality vertex from consideration
209
210 //first check if a muon is re-used in multiple vertices
211
212 if (workVerticesContainer.size() < 2) {
213 ATH_MSG_VERBOSE("Only one MuSA vertex, no need to select best!");
214 return StatusCode::SUCCESS;
215 }
216
217 //sort all the vertices by chi2 starting with the lowest
218 std::sort(workVerticesContainer.begin(), workVerticesContainer.end(), [](const MuSAVtxFitterTool::WrkVrt& v1, const MuSAVtxFitterTool::WrkVrt& v2) { return v1.chi2 < v2.chi2; });
219 ATH_MSG_VERBOSE("Sorted vertices by chi2! Lowest chi2: " << workVerticesContainer.front().chi2 << " Highest chi2: " << workVerticesContainer.back().chi2);
220
221 //starting with the lowest chi2 vertex, check its muons and remove any other vertices that use the same muons
222 for (unsigned int i = 0; i < workVerticesContainer.size(); i++) {
223 auto& vertex = workVerticesContainer.at(i);
224 if (!vertex.isGood) {
225 continue; // skip vertices we've already marked as bad
226 }
227 for (unsigned int j = i+1; j < workVerticesContainer.size(); j++) {
228 auto& otherVertex = workVerticesContainer.at(j);
229 for (const auto& muon : vertex.muonCandidates) {
230 if (std::find(otherVertex.muonCandidates.begin(), otherVertex.muonCandidates.end(), muon) != otherVertex.muonCandidates.end()) {
231 //remove the other vertex
232 ATH_MSG_VERBOSE("Found a vertex re-using a muon with chi2: " << otherVertex.chi2 << " removing from consideration!");
233 otherVertex.isGood = false;
234 }
235 }
236 }
237 }
238
239 //remove any vertices marked as bad
240 workVerticesContainer.erase(std::remove_if(workVerticesContainer.begin(), workVerticesContainer.end(),
241 [&](const MuSAVtxFitterTool::WrkVrt& vtx) { return !vtx.isGood; }), workVerticesContainer.end());
242
243 return StatusCode::SUCCESS;
244
245
246}
247} // namespace Rec
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
std::pair< std::vector< unsigned int >, bool > res
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
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)
Class describing the Line to which the Perigee refers to.
float beamPosX() const
X coordinate of the beam spot position.
float beamPosZ() const
Z coordinate of the beam spot position.
float beamPosY() const
Y coordinate of the beam spot position.
const Trk::Perigee & perigeeParameters() const
Returns the Trk::MeasuredPerigee track parameters.
virtual double pt() const override final
The transverse momentum ( ) of the particle.
virtual double eta() const override final
The pseudorapidity ( ) of the particle.
Eigen::Matrix< double, 3, 1 > Vector3D
Gaudi Tools.
PropDirection
PropDirection, enum for direction of the propagation.
@ oppositeMomentum
@ theta
Definition ParamDefs.h:66
@ qOverP
perigee
Definition ParamDefs.h:67
@ phi
Definition ParamDefs.h:75
@ d0
Definition ParamDefs.h:63
@ z0
Definition ParamDefs.h:64
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.
DataModel_detail::iterator< DVL > remove_if(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end, Predicate pred)
Specialization of remove_if for DataVector/List.
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.
double chi2
VKalVrt fit covariance.
std::vector< double > cov
VKalVrt fit vertex 4-momentum.
long int charge
list of VKalVrt fit chi2 for each track
double minOpAng
list of track parameters wrt the reconstructed vertex
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,...