ATLAS Offline Software
Loading...
Searching...
No Matches
GNNVertexFitterTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
6#include "InDetTrackSystematicsTools/InDetTrackTruthOriginDefs.h" //InDet::ExclusiveOrigin
7#include "GeoPrimitives/GeoPrimitivesHelpers.h" //Amg::deltaR
8#include "AthLinks/ElementLink.h"
9#include "xAODTracking/TrackParticleContainer.h" //template param for ElementLink
10
13#include <cmath>
14
15
16namespace Rec {
17
18GNNVertexFitterTool::GNNVertexFitterTool(const std::string &type, const std::string &name,
19 const IInterface *parent)
20 : AthAlgTool(type, name, parent), m_vertexFitterTool("Trk::TrkVKalVrtFitter/VertexFitterTool", this),
21 m_deco_mass("mass"),
22 m_deco_pt("pt"),
23 m_deco_charge("charge"),
24 m_deco_vPos("vPos"),
25 m_deco_lxy("Lxy"),
26 m_deco_sig3D("significance3D"),
27 m_deco_deltaR("deltaR"),
28 m_deco_ntrk("ntrk"),
29 m_deco_lxyz("Lxyz"),
30 m_deco_eFrac("efracsv"),
31 m_deco_nHFTracks("nHFtrk"){
32 declareInterface<IGNNVertexFitterInterface>(this);
33 declareProperty("JetTrackLinks", m_trackLinksKey = "AntiKt4EMPFlow."+m_gnnModel+"_TrackLinks");
34 declareProperty("JetTrackOrigins", m_trackOriginsKey = "AntiKt4EMPFlow."+m_gnnModel+"_TrackOrigin");
35 declareProperty("JetVertexLinks", m_vertexLinksKey = "AntiKt4EMPFlow."+m_gnnModel+"_VertexIndex");
36 declareProperty("VertexFitterTool", m_vertexFitterTool, "Vertex fitting tool");
38}
39
40/* Destructor */
41GNNVertexFitterTool::~GNNVertexFitterTool() { ATH_MSG_DEBUG("GNNVertexFitterTool destructor called"); }
42
44
45 ATH_MSG_DEBUG("GNNVertexFitter Tool in initialize()");
46
47 // Initialize keys
48 ATH_CHECK(m_trackLinksKey.initialize());
49 ATH_CHECK(m_trackOriginsKey.initialize());
50 ATH_CHECK(m_vertexLinksKey.initialize());
51
52 std::string linkNameMod = "";
53 if (m_doInclusiveVertexing) linkNameMod = "Inclusive";
54 m_jetWriteDecorKeyVertexLink = m_jetCollection + "."+linkNameMod+"GNNVerticesLink";
56
57 // Retrieve tools
58 ATH_CHECK(m_vertexFitterTool.retrieve());
59 // Additional Info for Vertex Fit
60 ATH_CHECK(m_beamSpotKey.initialize());
61 ATH_CHECK(m_eventInfoKey.initialize());
62
63 return StatusCode::SUCCESS;
64}
65
66// Total Momentum of Jet
67TLorentzVector GNNVertexFitterTool::TotalMom(const std::vector<const xAOD::TrackParticle *> &selTrk) const {
68 TLorentzVector sum(0., 0., 0., 0.);
69 for (int i = 0; i < (int)selTrk.size(); ++i) {
70 if (selTrk[i] == NULL)
71 continue;
72 sum += selTrk[i]->p4();
73 }
74 return sum;
75}
76
77// Vertex to Vertex Distance
78double GNNVertexFitterTool::vrtVrtDist(const xAOD::Vertex &primVrt, const Amg::Vector3D &secVrt,
79 const std::vector<double> &secVrtErr, double &signif) const {
80 double distx = primVrt.x() - secVrt.x();
81 double disty = primVrt.y() - secVrt.y();
82 double distz = primVrt.z() - secVrt.z();
83
84 AmgSymMatrix(3) primCovMtx = primVrt.covariancePosition(); // Create
85 primCovMtx(0, 0) += secVrtErr[0];
86 primCovMtx(0, 1) += secVrtErr[1];
87 primCovMtx(1, 0) += secVrtErr[1];
88 primCovMtx(1, 1) += secVrtErr[2];
89 primCovMtx(0, 2) += secVrtErr[3];
90 primCovMtx(2, 0) += secVrtErr[3];
91 primCovMtx(1, 2) += secVrtErr[4];
92 primCovMtx(2, 1) += secVrtErr[4];
93 primCovMtx(2, 2) += secVrtErr[5];
94
95 AmgSymMatrix(3) wgtMtx = primCovMtx.inverse();
96
97 signif = distx * wgtMtx(0, 0) * distx + disty * wgtMtx(1, 1) * disty + distz * wgtMtx(2, 2) * distz +
98 2. * distx * wgtMtx(0, 1) * disty + 2. * distx * wgtMtx(0, 2) * distz + 2. * disty * wgtMtx(1, 2) * distz;
99 signif = std::sqrt(std::abs(signif));
100 if (signif != signif)
101 signif = 0.;
102 return std::sqrt(distx * distx + disty * disty + distz * distz);
103}
104
105// Perform Vertex fit using the jet decorations of the GNN
107 xAOD::VertexContainer *outVertexContainer,
108 const xAOD::Vertex &primVrt, const EventContext &ctx) const {
109
110 using TLC = std::vector<ElementLink<DataVector<xAOD::TrackParticle_v1>>>;
112
113 // Read Decor Handle for Track links and Vertex links
118 jetWriteDecorHandleVertexLink(m_jetWriteDecorKeyVertexLink, ctx);
119
120 // Loop over the jets
121 for (const xAOD::Jet* jet : *inJetContainer) {
122
123 // Retrieve the Vertex and Track Collections
124 auto vertexCollection = vertexLinksHandle(*jet);
125 auto trackCollection = trackLinksHandle(*jet);
126 auto trackOriginCollection = trackOriginsHandle(*jet);
127
128 using indexList = std::vector<int>;
129 using trackCountMap = std::map<char, std::set<TL>>;
130
131 indexList iList(vertexCollection.size());
132
133 trackCountMap allTracksMap; // All the vertices and the corresponding track links
134 trackCountMap heavyFlavourVertexMap; // All Heavy Flavour Tracks associated with a vertex
135 trackCountMap primaryVertexMap; // All Primary Tracks associated with a vertex
136 trackCountMap fittingMap; // Map filled with vertices to be fitted
137
138 fittingMap.clear();
139 heavyFlavourVertexMap.clear();
140 primaryVertexMap.clear();
141 allTracksMap.clear();
142
143 // Fill the map of predicted vertex indices -> tracks
144 for (int index=0; index < int(vertexCollection.size()); index++){
145
146 auto vertex = vertexCollection[index];
147 auto trackOrigin = trackOriginCollection[index];
148 auto trackLink = trackCollection[index];
149
150 allTracksMap[vertex].insert(trackLink);
151
152 // Add heavy flavour tracks associated to each vertex to the map
153 if (InDet::ExclusiveOrigin::FromB == trackOrigin || InDet::ExclusiveOrigin::FromBC == trackOrigin ||
154 InDet::ExclusiveOrigin::FromC == trackOrigin) {
155 heavyFlavourVertexMap[vertex].insert(trackLink);
156 }
157 // Add primary tracks associated to each vertex to the map
158 else if (InDet::ExclusiveOrigin::Primary == trackOrigin) {
159 primaryVertexMap[vertex].insert(trackLink);
160 }
161 };
162
163 // determine the vertex with the largest number of primary tracks
164 auto pvCandidate = std::max_element(
165 primaryVertexMap.begin(), primaryVertexMap.end(),
166 [](const auto& a, const auto& b) { return a.second.size() < b.second.size(); }
167 )->first;
168
169 // filter the vertices according to the configurable options
170 for (const auto &vertexTrackPair : allTracksMap) {
171 // unpack vertexTrackPair
172 const auto &[vertex, trackLinks] = vertexTrackPair;
173 // remove the primary vertex unless requested
175 if(vertex == pvCandidate) {
176 continue;
177 }
178 }
179 // remove vertices which do not contain any heavy flavour tracks if requested
181 if (heavyFlavourVertexMap.find(vertex) == heavyFlavourVertexMap.end()) {
182 continue;
183 }
184 }
185 // remove vertices which contain insufficient heavy flavour tracks if requested
186 if(m_HFRatioThres > 0) {
187 if (heavyFlavourVertexMap.find(vertex) != heavyFlavourVertexMap.end() &&
188 (static_cast<float>(heavyFlavourVertexMap[vertex].size()) / trackLinks.size()) < m_HFRatioThres) {
189 continue;
190 }
191 }
192 // now we've passed all vertex filtering, add the vertex to the fittingMap
193 fittingMap.insert(std::pair<char, std::set<TL>>(vertex, trackLinks));
194 }
195 // If inclusive vertexing is requested, merge the vertices
197 // Define the key for the merged vertex
198 const char mergedVertexKey = 0;
199 std::set<TL> mergedSet;
200 std::set<TL> combinedHFTracks;
201
202 // Iterate over the fittingMap to merge all sets into a single set
203 for (const auto &vertexTrackPair : fittingMap) {
204 // Insert all elements from trackLinks into mergedSet
205 mergedSet.insert(vertexTrackPair.second.begin(), vertexTrackPair.second.end());
206 combinedHFTracks.insert(heavyFlavourVertexMap[vertexTrackPair.first].begin(), heavyFlavourVertexMap[vertexTrackPair.first].end());
207 }
208
209 // Clear fittingMap and insert the merged set under the mergedVertexKey
210 fittingMap.clear();
211 fittingMap[mergedVertexKey] = std::move(mergedSet);
212 heavyFlavourVertexMap.clear();
213 heavyFlavourVertexMap[mergedVertexKey] = std::move(combinedHFTracks);
214 }
215
216 // Working xAOD
217 std::unique_ptr<workVectorArrxAOD> xAODwrk (new workVectorArrxAOD());
219
220 // Beam Conditions
221 xAODwrk->beamX = beamSpotHandle->beamPos().x();
222 xAODwrk->beamY = beamSpotHandle->beamPos().y();
223 xAODwrk->beamZ = beamSpotHandle->beamPos().z();
224 xAODwrk->tanBeamTiltX = tan(beamSpotHandle->beamTilt(0));
225 xAODwrk->tanBeamTiltY = tan(beamSpotHandle->beamTilt(1));
226
227 std::unique_ptr<std::vector<WrkVrt>> wrkVrtSet = std::make_unique<std::vector<WrkVrt>>();
228 WrkVrt newvrt;
229 newvrt.Good = true;
230 std::unique_ptr<Trk::IVKalState> state = m_vertexFitterTool->makeState();
231 std::vector<const xAOD::NeutralParticle *> neutralPartDummy(0);
232 Amg::Vector3D IniVrt(0., 0., 0.);
233
234 for (const auto &pair : fittingMap) {
235 // Need at least 2 tracks to perform a fit
236 if (pair.second.size() >= 2) {
237 int NTRKS = pair.second.size();
238 std::vector<double> InpMass(NTRKS, m_massPi);
239 m_vertexFitterTool->setMassInputParticles(InpMass, *state);
240
241 xAODwrk->listSelTracks.clear();
242
243 for (auto TrackLink : pair.second) {
244 xAODwrk->listSelTracks.push_back(*TrackLink);
245 }
246
247 Amg::Vector3D FitVertex, vDist;
248 TLorentzVector jetDir(jet->p4()); // Jet Direction
249
250 // Get Estimate
251 StatusCode sc = m_vertexFitterTool->VKalVrtFitFast(xAODwrk->listSelTracks, FitVertex, *state);
252
253 if (sc.isFailure() || FitVertex.perp() > m_maxLxy) { /* No initial estimation */
254 IniVrt = primVrt.position();
256 IniVrt.setZero();
257 } else {
258 vDist = FitVertex - primVrt.position();
259 double JetVrtDir = jetDir.Px() * vDist.x() + jetDir.Py() * vDist.y() + jetDir.Pz() * vDist.z();
261 JetVrtDir = fabs(JetVrtDir); /* Always positive when primary vertex is seeked for*/
262 if (JetVrtDir > 0.)
263 IniVrt = FitVertex; /* Good initial estimation */
264 else
265 IniVrt = primVrt.position();
266 }
267
268 m_vertexFitterTool->setApproximateVertex(IniVrt.x(), IniVrt.y(), IniVrt.z(), *state);
269
270 // Perform the Vertex Fit
271 sc = (m_vertexFitterTool->VKalVrtFit(xAODwrk->listSelTracks, neutralPartDummy, newvrt.vertex, newvrt.vertexMom,
272 newvrt.vertexCharge, newvrt.vertexCov, newvrt.chi2PerTrk, newvrt.trkAtVrt,
273 newvrt.chi2, *state, false));
274 if (sc.isFailure())
275 continue;
276
277 // Chi2 Cut
278 auto NDOF = 2 * (newvrt.trkAtVrt.size()) - 3.0; // From VrtSecInclusive
279
280 if (newvrt.chi2 / NDOF >= m_maxChi2)
281 continue;
282 ATH_MSG_DEBUG("Found IniVertex=" << newvrt.vertex[0] << ", " << newvrt.vertex[1] << ", " << newvrt.vertex[2]
283 << " trks " << newvrt.trkAtVrt.size());
284
285 Amg::Vector3D vDir = newvrt.vertex - primVrt.position(); // Vertex Dirction in relation to Primary
286
287 Amg::Vector3D jetVrtDir(jet->p4().Px(), jet->p4().Py(), jet->p4().Pz());
288
289 double vPos =
290 (vDir.x() * newvrt.vertexMom.Px() + vDir.y() * newvrt.vertexMom.Py() + vDir.z() * newvrt.vertexMom.Pz()) /
291 newvrt.vertexMom.Rho();
292
293 double Lxy = sqrt(vDir[0] * vDir[0] + vDir[1] * vDir[1]);
294 double Lxyz = sqrt(vDir[0] * vDir[0] + vDir[1] * vDir[1] + vDir[2] * vDir[2]);
295 ATH_MSG_DEBUG("Lxyz " << Lxyz);
296
297 double drJPVSV = Amg::deltaR(jetVrtDir, vDir); // DeltaR
298
299 int ntrk = newvrt.trkAtVrt.size(); // # Tracks in Vertex
300
301 TLorentzVector MomentumVtx = TotalMom(xAODwrk->listSelTracks);
302
303 double eRatio = MomentumVtx.E() / jet->p4().E();
304 double signif3D;
305 [[maybe_unused]] double distToPV = vrtVrtDist(primVrt, newvrt.vertex, newvrt.vertexCov, signif3D);
306
307 // apply quality cuts
308 if (m_applyCuts) {
309 // cut on minimum number of tracks in the vertex
310 if (ntrk < m_minNTrack)
311 continue;
312 // cut on 3D significance
313 if (signif3D < m_minSig3D)
314 continue;
315 // cut on minumum transverse displacement from the PV
316 if (Lxy < m_minLxy )
317 continue;
318 }
319
320 // Register Container
321 auto* GNNvertex = outVertexContainer->emplace_back(new xAOD::Vertex);
322
323 for (const auto *trk : xAODwrk->listSelTracks) {
325 *(dynamic_cast<const xAOD::TrackParticleContainer *>(trk->container())),
326 trk->index());
327 GNNvertex->addTrackAtVertex(link_trk, 1.);
328 }
329
330 //Add Vertex Info into Container
332 GNNvertex->setPosition(newvrt.vertex);
333 GNNvertex->setFitQuality(newvrt.chi2, NDOF);
334
335 m_deco_mass(*GNNvertex) = newvrt.vertexMom.M();
336 m_deco_pt(*GNNvertex) = newvrt.vertexMom.Perp();
337 m_deco_charge(*GNNvertex) = newvrt.vertexCharge;
338 m_deco_vPos(*GNNvertex) = vPos;
339 m_deco_lxy(*GNNvertex) = Lxy;
340 m_deco_lxyz(*GNNvertex) = Lxyz;
341 m_deco_sig3D(*GNNvertex) = signif3D;
342 m_deco_ntrk(*GNNvertex) = ntrk;
343 m_deco_deltaR(*GNNvertex) = drJPVSV;
344 m_deco_eFrac(*GNNvertex) = eRatio;
345 m_deco_nHFTracks(*GNNvertex) = heavyFlavourVertexMap[pair.first].size();
346
348 linkVertex.setElement(GNNvertex);
349 linkVertex.setStorableObject(*outVertexContainer);
350 jetWriteDecorHandleVertexLink(*jet).push_back(linkVertex);
351
352 } // end of 2 Track requirement
353 }
354 } // end loop over jets
355 return StatusCode::SUCCESS;
356} // end performVertexFit
357
358StatusCode GNNVertexFitterTool::finalize() { return StatusCode::SUCCESS; }
359
360} // namespace Rec
361
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_DEBUG(x)
#define AmgSymMatrix(dim)
static Double_t a
static Double_t sc
ElementLink< xAOD::TrackParticleContainer > TrackLink
Handle class for reading a decoration on an object.
Handle class for adding a decoration to an object.
AthAlgTool(const std::string &type, const std::string &name, const IInterface *parent)
Constructor with parameters:
Gaudi::Details::PropertyBase & declareProperty(Gaudi::Property< T, V, H > &t)
value_type emplace_back(value_type pElem)
Add an element to the end of the collection.
SG::WriteDecorHandleKey< xAOD::JetContainer > m_jetWriteDecorKeyVertexLink
ToolHandle< Trk::TrkVKalVrtFitter > m_vertexFitterTool
SG::ReadCondHandleKey< InDet::BeamSpotData > m_beamSpotKey
SG::ReadHandleKey< xAOD::EventInfo > m_eventInfoKey
SG::AuxElement::Decorator< float > m_deco_deltaR
BooleanProperty m_removeNonHFVertices
SG::AuxElement::Decorator< float > m_deco_eFrac
BooleanProperty m_doInclusiveVertexing
SG::ReadDecorHandleKey< xAOD::JetContainer > m_trackOriginsKey
SG::AuxElement::Decorator< float > m_deco_pt
TLorentzVector TotalMom(const std::vector< const xAOD::TrackParticle * > &selTrk) const
SG::AuxElement::Decorator< float > m_deco_charge
SG::AuxElement::Decorator< float > m_deco_lxyz
BooleanProperty m_includePrimaryVertex
SG::ReadDecorHandleKey< xAOD::JetContainer > m_trackLinksKey
GNNVertexFitterTool(const std::string &type, const std::string &name, const IInterface *parent)
SG::AuxElement::Decorator< float > m_deco_vPos
SG::AuxElement::Decorator< float > m_deco_nHFTracks
SG::AuxElement::Decorator< float > m_deco_ntrk
double vrtVrtDist(const xAOD::Vertex &primVrt, const Amg::Vector3D &secVrt, const std::vector< double > &vrtErr, double &signif) const
SG::AuxElement::Decorator< float > m_deco_sig3D
SG::ReadDecorHandleKey< xAOD::JetContainer > m_vertexLinksKey
virtual StatusCode fitAllVertices(const xAOD::JetContainer *, xAOD::VertexContainer *, const xAOD::Vertex &primaryVertex, const EventContext &) const
SG::AuxElement::Decorator< float > m_deco_mass
SG::AuxElement::Decorator< float > m_deco_lxy
Handle class for reading a decoration on an object.
Handle class for adding a decoration to an object.
STL class.
float z() const
Returns the z position.
void setVertexType(VxType::VertexType vType)
Set the type of the vertex.
float y() const
Returns the y position.
const Amg::Vector3D & position() const
Returns the 3-pos.
float x() const
Returns the x position.
Eigen::Matrix< double, 3, 1 > Vector3D
double deltaR(const Amg::Vector3D &v1, const Amg::Vector3D &v2)
constexpr double chargedPionMassInMeV
the mass of the charged pion (in MeV)
Gaudi Tools.
Definition index.py:1
@ SecVtx
Secondary vertex.
Jet_v1 Jet
Definition of the current "jet version".
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".
JetContainer_v1 JetContainer
Definition of the current "jet container version".