ATLAS Offline Software
Loading...
Searching...
No Matches
SiSPGNNTrackMaker.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5#include <memory>
6#include <fstream>
7
8#include "SiSPGNNTrackMaker.h"
9
12
14 const std::string& name, ISvcLocator* pSvcLocator)
15 : AthReentrantAlgorithm(name, pSvcLocator)
16{
17
18}
19
21{
22 ATH_CHECK(m_SpacePointsPixelKey.initialize());
23 ATH_CHECK(m_SpacePointsSCTKey.initialize());
25 ATH_CHECK(m_ClusterPixelKey.initialize());
26 ATH_CHECK(m_ClusterStripKey.initialize());
27
28 ATH_CHECK(m_outputTracksKey.initialize());
29
30 ATH_CHECK(m_trackFitter.retrieve());
31 ATH_CHECK(m_seedFitter.retrieve());
32 ATH_CHECK(m_trackSummaryTool.retrieve());
33
34 if (!m_gnnTrackFinder.empty() && !m_gnnTrackReader.empty()) {
35 ATH_MSG_ERROR("Use either track finder or track reader, not both.");
36 return StatusCode::FAILURE;
37 }
38
39 if (!m_gnnTrackFinder.empty()) {
40 ATH_MSG_INFO("Use GNN Track Finder");
41 ATH_CHECK(m_gnnTrackFinder.retrieve());
42 }
43 if (!m_gnnTrackReader.empty()) {
44 ATH_MSG_INFO("Use GNN Track Reader");
45 ATH_CHECK(m_gnnTrackReader.retrieve());
46 }
48 ATH_MSG_INFO("Use input clusters");
49 }
50 // retrieve eta dependent cut svc
52
53 ATH_MSG_INFO("Applying the following cuts during GNN-based track reconstruction: ");
54 ATH_MSG_INFO("Min pT: " << m_pTmin);
55 ATH_MSG_INFO("Max eta: " << m_etamax);
56 ATH_MSG_INFO("Min number of clusters: " << m_minClusters);
57 ATH_MSG_INFO("Min number of pixel clusters: " << m_minPixelClusters);
58 ATH_MSG_INFO("Min number of strip clusters: " << m_minStripClusters);
59
60
62 ATH_MSG_INFO("Applying the following eta dependant track cuts after GNN-based track reconstruction: ");
63 std::vector <double> etaBins, minPT, maxz0, maxd0;
64 std::vector <int> minClusters, minPixelHits, maxHoles;
65 m_etaDependentCutsSvc->getValue(InDet::CutName::etaBins, etaBins);
66 ATH_MSG_INFO("Eta bins: " << etaBins );
67 m_etaDependentCutsSvc->getValue(InDet::CutName::minClusters, minClusters);
68 ATH_MSG_INFO("Min Si Hits: " << minClusters );
69 m_etaDependentCutsSvc->getValue(InDet::CutName::minPixelHits, minPixelHits);
70 ATH_MSG_INFO("Min pixel hits: " << minPixelHits );
71 m_etaDependentCutsSvc->getValue(InDet::CutName::maxHoles, maxHoles);
72 ATH_MSG_INFO("Max holes: " << maxHoles);
73 m_etaDependentCutsSvc->getValue(InDet::CutName::minPT, minPT);
74 ATH_MSG_INFO("Min pT: " << minPT);
75 m_etaDependentCutsSvc->getValue(InDet::CutName::maxZImpact, maxz0);
76 ATH_MSG_INFO("Max z0: " << maxz0);
77 m_etaDependentCutsSvc->getValue(InDet::CutName::maxPrimaryImpact, maxd0);
78 ATH_MSG_INFO("Max d0: " << maxd0);
79 }
80
81 return StatusCode::SUCCESS;
82}
83
84StatusCode InDet::SiSPGNNTrackMaker::execute(const EventContext& ctx) const
85{
87 ATH_CHECK(outputTracks.record(std::make_unique<TrackCollection>()));
88
89 // get event info
90 uint32_t runNumber = ctx.eventID().run_number();
91 uint32_t eventNumber = ctx.eventID().event_number();
92
93 // get all space points from event
94 std::vector<const Trk::SpacePoint*> spacePoints = getSpacePointsInEvent(ctx, eventNumber);
95
96 // get all clusters from event
97 std::vector<const Trk::PrepRawData*> allClusters = getClustersInEvent(ctx, eventNumber);
98
99 // get tracks from GNN chain
100 std::vector<std::vector<uint32_t> > TT;
101 std::vector<std::vector<uint32_t> > clusterTracks;
102 if (m_gnnTrackFinder.isSet()) {
103 ATH_CHECK(m_gnnTrackFinder->getTracks(spacePoints, TT));
104 } else if (m_gnnTrackReader.isSet()) {
105 // if track candidates are built from cluster, get both clusters and SPs
107 m_gnnTrackReader->getTracks(runNumber, eventNumber, clusterTracks, TT) :
108 m_gnnTrackReader->getTracks(runNumber, eventNumber, TT);
109 } else {
110 ATH_MSG_ERROR("Both GNNTrackFinder and GNNTrackReader are not set");
111 return StatusCode::FAILURE;
112 }
113
114 ATH_MSG_DEBUG("Event " << eventNumber << " obtained " << TT.size() << " Tracks");
115
116 // loop over all track candidates
117 // and perform track fitting for each.
118 int trackCounter = -1;
119 std::vector<int> status_codes;
120 // track processing loop
121 for (auto& trackIndices : TT) {
122 // For each track candidate:
123 trackCounter++;
124
125 // 1. Sort space points by distance from origin
126 std::vector<const Trk::SpacePoint*> trackCandidate = getSpacePoints(trackIndices, spacePoints);
127
128 // 2. Get associated clusters
129 // if track candidates are built from cluster, get both clusters and SPs
130 std::vector<const Trk::PrepRawData*> clusters = m_areInputClusters ?
131 getClusters (clusterTracks, allClusters, trackCounter)
132 : spacePointsToClusters(trackCandidate) ;
133
134 // 3. Perform track fitting:
135 // - Initial conformal mapping
136 // - First chi2 fit without outlier removal
137 // - Second chi2 fit with improved initial estimate of track parameters and with outlier removal
138 // 4. Apply quality cuts (pT, eta)
139 // 5. Compute track summary
140 // 6. Store track if it passes all criteria
141 auto [fitSuccess, passTrackCut, track] = doFitAndCut(ctx, trackCandidate, clusters, trackCounter);
142
143 if (not fitSuccess) {
144 continue;
145 }
146
147 if (passTrackCut <= 0) {
148 outputTracks->push_back(track.release());
149 }
150
151 status_codes.push_back(passTrackCut);
152
153 }
154
155 if (m_doRecoTrackCuts) {
156 ATH_MSG_INFO("Event " << eventNumber << " has " << status_codes.size() << " tracks found, "
157 << std::count(status_codes.begin(), status_codes.end(), 0)
158 << " tracks remains after applying track cuts");
159 } else {
160 ATH_MSG_INFO("Event " << eventNumber << " has " << status_codes.size() << " tracks found, all tracks are kept");
161 }
162
163 return StatusCode::SUCCESS;
164}
165
166
167bool InDet::SiSPGNNTrackMaker::prefitCheck(unsigned int nPix, unsigned int nStrip, unsigned int nClusters, unsigned int nSpacePoints) const {
168 return nPix >= m_minPixelClusters && nStrip >= m_minStripClusters && nClusters >= m_minClusters && nSpacePoints >= 3;
169}
170
171
173{
174 const Trk::Perigee* origPerigee = track.perigeeParameters();
175 double pt = origPerigee->pT();
176
177 double eta = std::abs(origPerigee->eta());
178
179 double d0 = std::abs(origPerigee->parameters()[Trk::d0]);
180
181 double z0 = std::abs(origPerigee->parameters()[Trk::z0]);
182
183 int nHolesOnTrack = track.trackSummary()->get(Trk::numberOfPixelHoles) +
184 track.trackSummary()->get(Trk::numberOfSCTHoles);
185
186 int nPixels = track.trackSummary()->get(Trk::numberOfPixelHits);
187 int nStrips = track.trackSummary()->get(Trk::numberOfSCTHits);
188 int nClusters = nPixels + nStrips;
189
190 ATH_MSG_DEBUG("track params: " << pt << " " << eta << " " << d0 << " " << z0
191 << " " << nClusters << nStrips
192 << " " << nPixels);
193
194 // min Si hits
195 if (nClusters < m_etaDependentCutsSvc->getMinSiHitsAtEta(eta))
196 return 1;
197
198 // min pixel hits
199 if (nPixels < m_etaDependentCutsSvc->getMinPixelHitsAtEta(eta))
200 return 2;
201
202 // min pT, default 400
203 if (pt < m_etaDependentCutsSvc->getMinPtAtEta(eta))
204 return 3;
205
206 // max z0
207 if (z0 > m_etaDependentCutsSvc->getMaxZImpactAtEta(eta))
208 return 4;
209
210 // max d0
211 if (d0 > m_etaDependentCutsSvc->getMaxPrimaryImpactAtEta(eta))
212 return 5;
213
214 // max holes
215 if (nHolesOnTrack > m_etaDependentCutsSvc->getMaxSiHolesAtEta(eta))
216 return 6;
217
218 return 0;
219}
220
221std::vector<const Trk::SpacePoint*> InDet::SiSPGNNTrackMaker::getSpacePointsInEvent (
222 const EventContext& ctx,
223 int eventNumber
224) const {
225 std::vector<const Trk::SpacePoint*> spacePoints;
226
227 int npixsp(0), nstrip(0), n_overlap(0);
229 spacePoints.push_back(sp);
230 npixsp++;
231 }
233 spacePoints.push_back(sp);
234 nstrip++;
235 }
237 spacePoints.push_back(sp);
238 n_overlap++;
239 }
240 ATH_MSG_DEBUG("Event " << eventNumber << " has " << npixsp
241 << " pixel space points, " << nstrip
242 << " strips space points" << n_overlap
243 << " overlapping spacepoints" << spacePoints.size()
244 << " space points");
245
246 return spacePoints;
247}
248
249std::vector<const Trk::SpacePoint*> InDet::SiSPGNNTrackMaker::getSpacePointsInEvent(
250 const EventContext& ctx,
252) const {
253 std::vector<const Trk::SpacePoint*> spacePoints;
254 if (not containerKey.empty()){
255
257
258 if (container.isValid()){
259 // loop over spacepoint collection
260 auto spc = container->begin();
261 auto spce = container->end();
262 for(; spc != spce; ++spc){
263 const SpacePointCollection* spCollection = (*spc);
264 auto sp = spCollection->begin();
265 auto spe = spCollection->end();
266 for(; sp != spe; ++sp) {
267 const Trk::SpacePoint* spacePoint = (*sp);
268 spacePoints.push_back(spacePoint);
269 }
270 }
271 }
272 }
273
274 return spacePoints;
275}
276
277std::vector<const Trk::SpacePoint*> InDet::SiSPGNNTrackMaker::getSpacePointsInEvent(
278 const EventContext& ctx,
280) const {
281
282 std::vector<const Trk::SpacePoint*> spacePoints;
283
284 if (not containerKey.empty()){
285
286 SG::ReadHandle<SpacePointOverlapCollection> collection{containerKey, ctx};
287
288 if (collection.isValid()){
289 for (const Trk::SpacePoint *sp : *collection) {
290 spacePoints.push_back(sp);
291 }
292 }
293 }
294
295 return spacePoints;
296}
297
298std::vector<const Trk::PrepRawData*> InDet::SiSPGNNTrackMaker::getClustersInEvent (
299 const EventContext& ctx,
300 int eventNumber
301) const {
302 std::vector<const Trk::PrepRawData*> allClusters;
303 if (m_areInputClusters) {
305 ctx);
307 ctx);
308
309 if (!pixcontainer.isValid()) {
310 ATH_MSG_ERROR("Pixel container invalid, returning");
311 return allClusters;
312 }
313
314 if (!strip_container.isValid()) {
315 ATH_MSG_ERROR("Strip container invalid, returning");
316 return allClusters;
317 }
318
319 auto pixcollection = pixcontainer->begin();
320 auto pixcollectionEnd = pixcontainer->end();
321 for (; pixcollection != pixcollectionEnd; ++pixcollection) {
322 if ((*pixcollection)->empty()) {
323 ATH_MSG_WARNING("Empty pixel cluster collection encountered");
324 continue;
325 }
326 auto const* clusterCollection = (*pixcollection);
327 auto thisCluster = clusterCollection->begin();
328 auto clusterEnd = clusterCollection->end();
329 for (; thisCluster != clusterEnd; ++thisCluster) {
330 const PixelCluster* cl = (*thisCluster);
331 allClusters.push_back(cl);
332 }
333 }
334
335 auto strip_collection = strip_container->begin();
336 auto strip_collectionEnd = strip_container->end();
337 for (; strip_collection != strip_collectionEnd; ++strip_collection) {
338 if ((*strip_collection)->empty()) {
339 ATH_MSG_WARNING("Empty strip cluster collection encountered");
340 continue;
341 }
342 auto const* clusterCollection = (*strip_collection);
343 auto thisCluster = clusterCollection->begin();
344 auto clusterEnd = clusterCollection->end();
345 for (; thisCluster != clusterEnd; ++thisCluster) {
346 const SCT_Cluster* cl = (*thisCluster);
347 allClusters.push_back(cl);
348 }
349 }
350
351 ATH_MSG_DEBUG("Event " << eventNumber << " has " << allClusters.size()
352 << " clusters");
353 }
354 return allClusters;
355}
356
357
358std::vector<const Trk::SpacePoint*> InDet::SiSPGNNTrackMaker::getSpacePoints (
359 const std::vector<uint32_t>& trackIndices,
360 const std::vector<const Trk::SpacePoint*>& allSpacePoints
361) const {
362
363 std::vector<const Trk::SpacePoint*> trackCandidate;
364 trackCandidate.reserve(trackIndices.size());
365
366 std::vector<std::pair<double, const Trk::SpacePoint*> > distanceSortedSPs;
367
368 // get track space points
369 // sort SPs in track by distance from origin
370 for (auto& id : trackIndices) {
372 if (id > allSpacePoints.size()) {
373 ATH_MSG_WARNING("SpacePoint index "
374 << id << " out of range: " << allSpacePoints.size());
375 continue;
376 }
377
378 const Trk::SpacePoint* sp = allSpacePoints[id];
379
380 // store distance - hit paire
381 if (sp != nullptr) {
382 distanceSortedSPs.push_back(
383 std::make_pair(
384 pow(sp->globalPosition().x(), 2) + pow(sp->globalPosition().y(), 2),
385 sp
386 )
387 );
388 }
389 }
390
391 // sort by distance
392 std::sort(distanceSortedSPs.begin(), distanceSortedSPs.end());
393
394 // add SP to trk candidate in the same order
395 for (size_t i = 0; i < distanceSortedSPs.size(); i++) {
396 trackCandidate.push_back(distanceSortedSPs[i].second);
397 }
398
399 return trackCandidate;
400}
401
402std::vector<const Trk::PrepRawData*> InDet::SiSPGNNTrackMaker :: spacePointsToClusters (
403 const std::vector<const Trk::SpacePoint*>& spacePoints
404) const {
405 std::vector<const Trk::PrepRawData*> clusters;
406 for (const Trk::SpacePoint* sp : spacePoints) {
407 clusters.push_back(sp->clusterList().first);
408 if (sp->clusterList().second != nullptr) {
409 clusters.push_back(sp->clusterList().second);
410 }
411 }
412 return clusters;
413}
414
415std::vector<const Trk::PrepRawData*> InDet::SiSPGNNTrackMaker :: getClusters (
416 const std::vector<std::vector<uint32_t>>& clusterTracks,
417 const std::vector<const Trk::PrepRawData*>& allClusters,
418 int trackNumber
419) const {
420 std::vector<uint32_t> clusterIndices = clusterTracks[trackNumber];
421 std::vector<const Trk::PrepRawData*> clusters;
422 clusters.reserve(clusterIndices.size());
423 for (uint32_t id : clusterIndices) {
424 if (id > allClusters.size()) {
425 ATH_MSG_ERROR("Cluster index out of range");
426 continue;
427 }
428 clusters.push_back(allClusters[id]);
429 }
430 return clusters;
431}
432
465
466std::tuple<bool, int, std::unique_ptr<Trk::Track>> InDet::SiSPGNNTrackMaker::doFitAndCut (
467 const EventContext& ctx,
468 std::vector<const Trk::SpacePoint*>& spacePoints,
469 std::vector<const Trk::PrepRawData*>& clusters,
470 int& trackCounter
471 ) const
472 {
473 // get cluster list
474 int nPIX(0), nStrip(0);
475
476 for (const Trk::PrepRawData* cl : clusters) {
477 if (cl->type(Trk::PrepRawDataType::PixelCluster)) nPIX++;
478 if (cl->type(Trk::PrepRawDataType::SCT_Cluster)) nStrip++;
479 }
480
481 ATH_MSG_DEBUG("Track " << trackCounter << " has " << spacePoints.size()
482 << " space points, " << clusters.size()
483 << " clusters, " << nPIX << " pixel clusters, "
484 << nStrip << " strip clusters");
485
486 // check hit counts
487 if (not prefitCheck(nPIX, nStrip, clusters.size(), spacePoints.size())) {
488 ATH_MSG_DEBUG("Track " << trackCounter << " does not pass prefit cuts, skipping");
489 return std::make_tuple(false, 999, nullptr);
490 }
491
492 // conformal mapping for track parameters
493 auto trkParameters = m_seedFitter->fit(spacePoints);
494 if (trkParameters == nullptr) {
495 ATH_MSG_DEBUG("Conformal mapping failed");
496 return std::make_tuple(false, 999, nullptr);
497 }
498
499 std::unique_ptr<Trk::Track> track = fitTrack(ctx, clusters, *trkParameters, trackCounter);
500
501 if (track == nullptr || track->perigeeParameters() == nullptr) {
502 ATH_MSG_DEBUG("Track " << trackCounter
503 << " fails the chi2 fit, skipping");
504 return std::make_tuple(false, 999, nullptr);
505 }
506
507 // compute pT and skip if pT too low
508 if (track->perigeeParameters()->pT() < m_pTmin) {
509 ATH_MSG_DEBUG("Track " << trackCounter
510 << "with pt = " << track->perigeeParameters()->pT()
511 << " has pT too low, skipping track!");
512 return std::make_tuple(false, 999, nullptr);
513 }
514
515 // get rid of tracks with eta too large
516 if (std::abs(track->perigeeParameters()->eta()) > m_etamax) {
517 ATH_MSG_DEBUG("Track " << trackCounter << "with eta = "
518 << std::abs(track->perigeeParameters()->eta())
519 << " has eta too high, skipping track!");
520 return std::make_tuple(false, 999, nullptr);
521 }
522
523 // if track fit succeeds, eta and pT within range, compute track summary. This is quite expensive.
524 m_trackSummaryTool->computeAndReplaceTrackSummary(
525 *track, false /* DO NOT suppress hole search*/);
526
527 int passTrackCut = (m_doRecoTrackCuts) ? passEtaDepCuts(*track) : -1;
528
529 return std::make_tuple(true, passTrackCut, std::move(track));
530
531 }
532
533std::unique_ptr<Trk::Track> InDet::SiSPGNNTrackMaker::fitTrack (
534 const EventContext& ctx,
535 std::vector<const Trk::PrepRawData*> clusters,
536 const Trk::TrackParameters& initial_params,
537 int trackCounter
538 ) const
539 {
540
542 // first fit the track with local parameters and without outlier removal.
543 bool keepOnTrying = true;
544 std::unique_ptr<Trk::Track> track;
545 while (keepOnTrying) {
546 track = m_trackFitter->fit(ctx, clusters, initial_params, false, matEffects);
547 // any need to recover a failed fit ?
548 if (track == nullptr || track->perigeeParameters() == nullptr) {
549 clusters.pop_back();
550 if (clusters.size()<m_minClusters || !m_doRecoverFailedFits) {
551 keepOnTrying = false;
552 }
553 } else {
554 keepOnTrying = false;
555 }
556 }
557
558 if (track == nullptr || track->perigeeParameters() == nullptr) {
559 ATH_MSG_DEBUG("Track " << trackCounter
560 << " fails the first chi2 fit, skipping");
561 return track;
562 }
563
564 // reject track with pT too low, default 400 MeV
565 if (track->perigeeParameters()->pT() < m_pTmin) {
566 ATH_MSG_DEBUG("Track " << trackCounter
567 << " fails the first chi2 fit, skipping");
568 return nullptr;
569 }
570
571 // finally fit with outlier removal
572 Trk::Perigee origPerigee = *track->perigeeParameters();
573 keepOnTrying = true;
574 while (keepOnTrying) {
575 track = m_trackFitter->fit(ctx, clusters, origPerigee, true, matEffects);
576 // any need to recover a failed or bad fit ?
577 bool doRefit=false;
578 if (track == nullptr || track->trackSummary() == nullptr || track->outliersOnTrack()->size()>=3) {
579 doRefit=true;
580 }
581 if (doRefit && m_doRecoverFailedFits) {
582 clusters.pop_back();
583 if (clusters.size()<m_minClusters) {
584 keepOnTrying=false;
585 }
586 } else {
587 keepOnTrying=false;
588 }
589 }
590
591 if (!track) ATH_MSG_DEBUG("Track " << trackCounter
592 << " fails the second chi2 fit, skipping");
593
594 return track;
595
596 }
597
599// Overload of << operator MsgStream
601
602MsgStream& InDet::operator <<
603 (MsgStream& sl,const InDet::SiSPGNNTrackMaker& se)
604{
605 return se.dump(sl);
606}
607
609// Overload of << operator std::ostream
611std::ostream& InDet::operator <<
612 (std::ostream& sl,const InDet::SiSPGNNTrackMaker& se)
613{
614 return se.dump(sl);
615}
616
618// Dumps relevant information into the MsgStream
620
621MsgStream& InDet::SiSPGNNTrackMaker::dump( MsgStream& out ) const
622{
623 out<<std::endl;
624 if(msgLvl(MSG::DEBUG)) return dumpevent(out);
625 else return dumptools(out);
626}
627
629// Dumps conditions information into the MsgStream
631
632MsgStream& InDet::SiSPGNNTrackMaker::dumptools( MsgStream& out ) const
633{
634 out<<"| Location of output tracks | "
635 <<std::endl;
636 out<<"|----------------------------------------------------------------"
637 <<"----------------------------------------------------|"
638 <<std::endl;
639 return out;
640}
641
643// Dumps event information into the ostream
645
646MsgStream& InDet::SiSPGNNTrackMaker::dumpevent( MsgStream& out ) const
647{
648 return out;
649}
650
651std::ostream& InDet::SiSPGNNTrackMaker::dump( std::ostream& out ) const
652{
653 return out;
654}
Scalar eta() const
pseudorapidity method
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
static Double_t sp
constexpr int pow(int base, int exp) noexcept
bool msgLvl(const MSG::Level lvl) const
An algorithm that can be simultaneously executed in multiple threads.
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.
InDet::SiSPGNNTrackMaker is an algorithm that uses the GNN-based track finding tool to reconstruct tr...
ToolHandle< ISeedFitter > m_seedFitter
MsgStream & dump(MsgStream &out) const
std::unique_ptr< Trk::Track > fitTrack(const EventContext &ctx, std::vector< const Trk::PrepRawData * > clusters, const Trk::TrackParameters &initial_params, int trackCounter) const
UnsignedIntegerProperty m_minStripClusters
ToolHandle< IGNNTrackReaderTool > m_gnnTrackReader
SG::ReadHandleKey< InDet::PixelClusterContainer > m_ClusterPixelKey
BooleanProperty m_doRecoverFailedFits
std::vector< const Trk::SpacePoint * > getSpacePoints(const std::vector< uint32_t > &trackIndices, const std::vector< const Trk::SpacePoint * > &allSpacePoints) const
SG::ReadHandleKey< InDet::SCT_ClusterContainer > m_ClusterStripKey
UnsignedIntegerProperty m_minPixelClusters
BooleanProperty m_doRecoTrackCuts
BooleanProperty m_areInputClusters
int passEtaDepCuts(const Trk::Track &track) const
SiSPGNNTrackMaker(const std::string &name, ISvcLocator *pSvcLocator)
ServiceHandle< IInDetEtaDependentCutsSvc > m_etaDependentCutsSvc
std::vector< const Trk::PrepRawData * > getClustersInEvent(const EventContext &ctx, int eventNumber) const
SG::WriteHandleKey< TrackCollection > m_outputTracksKey
std::vector< const Trk::PrepRawData * > spacePointsToClusters(const std::vector< const Trk::SpacePoint * > &spacePoints) const
ToolHandle< Trk::ITrackFitter > m_trackFitter
Track Fitter.
SG::ReadHandleKey< SpacePointContainer > m_SpacePointsSCTKey
MsgStream & dumpevent(MsgStream &out) const
virtual StatusCode initialize() override
ToolHandle< IGNNTrackFinder > m_gnnTrackFinder
GNN-based track finding tool that produces track candidates.
bool prefitCheck(unsigned int nPix, unsigned int nStrip, unsigned int nClusters, unsigned int nSpacePoints) const
std::vector< const Trk::SpacePoint * > getSpacePointsInEvent(const EventContext &ctx, int eventNumber) const
UnsignedIntegerProperty m_minClusters
virtual StatusCode execute(const EventContext &ctx) const override
ToolHandle< Trk::IExtendedTrackSummaryTool > m_trackSummaryTool
std::vector< const Trk::PrepRawData * > getClusters(const std::vector< std::vector< uint32_t > > &clusterTracks, const std::vector< const Trk::PrepRawData * > &allClusters, int trackNumber) const
SG::ReadHandleKey< SpacePointOverlapCollection > m_SpacePointsOverlapKey
SG::ReadHandleKey< SpacePointContainer > m_SpacePointsPixelKey
std::tuple< bool, int, std::unique_ptr< Trk::Track > > doFitAndCut(const EventContext &ctx, std::vector< const Trk::SpacePoint * > &spacePoints, std::vector< const Trk::PrepRawData * > &clusters, int &trackCounter) const
Fits a track and applies quality cuts to determine if it should be kept.
MsgStream & dumptools(MsgStream &out) const
Property holding a SG store/key/clid from which a ReadHandle is made.
virtual bool isValid() override final
Can the handle be successfully dereferenced?
bool empty() const
Test if the key is blank.
StatusCode record(std::unique_ptr< T > data)
Record a const object to the store.
double eta() const
Access method for pseudorapidity - from momentum.
double pT() const
Access method for transverse momentum.
ParametersT< TrackParametersDim, Charged, PerigeeSurface > Perigee
@ d0
Definition ParamDefs.h:63
@ z0
Definition ParamDefs.h:64
ParticleHypothesis
Enumeration for Particle hypothesis respecting the interaction with material.
ParametersBase< TrackParametersDim, Charged > TrackParameters
@ numberOfPixelHits
number of pixel layers on track with absence of hits
@ numberOfSCTHoles
number of Holes in both sides of a SCT module
@ numberOfPixelHoles
number of pixels which have a ganged ambiguity.
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.