ATLAS Offline Software
Loading...
Searching...
No Matches
InDetDenseEnvAmbiTrackSelectionTool.h
Go to the documentation of this file.
1// -*- C++ -*-
2
3/*
4 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
5*/
6
8// InDetDenseEnvAmbiTrackSelectionTool.h, (c) ATLAS Detector software
10
11#ifndef INDETInDetDenseEnvAmbiTrackSelectionTool_H
12#define INDETInDetDenseEnvAmbiTrackSelectionTool_H
13
16
25#include "TrkTrack/Track.h" //for use in the struct lessTrkTrack implementation in this header
28
31
32#include "GaudiKernel/ToolHandle.h"
33
34#include <cmath> //for std::fabs in implementation of structs in this header
35#include <iostream> //for cout in dumpInfo
36#include <map>
37#include <mutex>
38#include <vector>
39
40class Identifier;
41class SiliconID;
42namespace Trk {
44 class PrepRawData;
45 class RIO_OnTrack;
46}
47
48namespace InDet
49{
50
51 //In here we refer to ROIs as region of interests corresponding to high pt B tracks.
60 class InDetDenseEnvAmbiTrackSelectionTool : public extends<AthAlgTool, Trk::IAmbiTrackSelectionTool>
61 {
62 public:
63 InDetDenseEnvAmbiTrackSelectionTool(const std::string&,const std::string&,const IInterface*);
64
67
69 virtual StatusCode initialize() override;
70 virtual StatusCode finalize() override;
71
81 virtual std::tuple<Trk::Track*,bool> getCleanedOutTrack(const Trk::Track *track,
82 const Trk::TrackScore score,
83 Trk::ClusterSplitProbabilityContainer &splitProbContainer,
84 Trk::PRDtoTrackMap &prd_to_track_map,
85 int trackId,
86 int subtrackId) const override;
87
88 private:
89
90 // possible classification of TSOS
91 enum TsosTypes {
92 // A measurement not yet used in any other track
94 // A measurement shared with another track that can be split
96 // A measurement shared with another track
98 // A hit that needs to be removed from the track
99 // because it is used on too many tracks already
101 // A hit that needs to be removed from the track
102 // because sharing it would invalidate an accepted track
104 // A hit that needs to be removed from the track
105 // reason other than the above two
107 // an outlier, to be copied in case
109 // other TSOS types to be copied in case
111 };
112
113
115 {
117 bool m_thisHasIBLHit ; //Track has a hit in the IBL
118 bool m_hasSharedIBLHit ; //Track has a shared hit in the IBL
119 bool m_hasSharedPixel ; //Track has a shared pixel hit (does not include IBL)
120 bool m_firstPixIsShared ; //The first pixel hit on track is shared (includes IBL)
121 bool m_passHadronicROI ; //Track is within hadronic ROI
122 bool m_passConversionSel ; //Track is compatible with a conversion
123 bool m_trkCouldBeAccepted ; //No rejected hits on track
124
125 int m_numPixelDeadSensor ; //Taken from track summary. Stored value not changed in this tool
126 int m_numSCTDeadSensor ; //Taken from track summary. Stored value not changed in this tool
127 int m_numPixelHits ; //Taken from track summary. Stored value not changed in this tool
128 int m_numPixelHoles ; //Taken from track summary. Stored value not changed in this tool
129 int m_numSCTHoles ; //Taken from track summary. Stored value not changed in this tool
130 int m_numSCTHits ; //Taken from track summary. Stored value not changed in this tool
131
132 // values counted in this tool
133 int m_numUnused ; //Number of unique hits
134 int m_numTRT_Unused ; //Number of unique TRT hits (NOT included within m_numUnused)
135 int m_numSCT_Unused ; //Number of unique SCT hits ( included within m_numUnused)
136 int m_numPseudo ; //Number of pseudo measurements on track
137
138 int m_numSplitSharedPix ; //Number of Pixel clusters comptaible with being split that are also shared
139 int m_numSplitSharedSCT ; //Number of SCT clusters comptaible with being split that are also shared
140 int m_numShared ; //Number of shared hits on track (does not include SplitShared)
141 int m_numSCT_Shared ; //Number of shared sct hits on track (included within m_numShared)
142 int m_numWeightedShared ; //Weighted number of shared hits on track used to test shared module cuts (pixels count for 2)
143
145 {
146 m_isPatternTrack = true;
147 m_thisHasIBLHit = false;
148 m_hasSharedIBLHit = false;
149 m_hasSharedPixel = false;
150 m_firstPixIsShared = false;
151 m_passHadronicROI = false;
152 m_passConversionSel = false;
154
157 m_numPixelHits = 0;
158 m_numPixelHoles = 0;
159 m_numSCTHoles = 0;
160 m_numSCTHits = 0;
161
162 m_numUnused = 0;
163 m_numSCT_Unused = 0;
164 m_numTRT_Unused = 0;
165 m_numPseudo = 0;
166
169 m_numShared = 0;
170 m_numSCT_Shared = 0;
172 }
173
174 // All SiHits including shared
179
180 // SplitShared are not counted as shared (no penality)
181 // since are compatible with originating from multiple particles
183 {
184 return (totalSiHits() - m_numShared);
185 };
186
187 // No counters, for pixel hits. Solve for them to avoid too many counters
189 {
190 // m_numTRT_Unused not included in m_numUnused
191 return (m_numUnused - m_numSCT_Unused);
192 };
193
194 // No counters, for pixel hits. Solve for them to avoid too many counters
196 {
197 return (m_numShared - m_numSCT_Shared);
198 };
199
204
205 void dumpInfo() const
206 {
207 std::cout << "isPatternTrack : " << m_isPatternTrack << std::endl;
208 std::cout << "thisHasIBLHit : " << m_thisHasIBLHit << std::endl;
209 std::cout << "hasSharedIBLHit : " << m_hasSharedIBLHit << std::endl;
210 std::cout << "hasSharedPixel : " << m_hasSharedPixel << std::endl;
211 std::cout << "firstPixIsShared : " << m_firstPixIsShared << std::endl;
212 std::cout << "passHadronicROI : " << m_passHadronicROI << std::endl;
213 std::cout << "passConversionSel : " << m_passConversionSel << std::endl;
214 std::cout << "trkCouldBeAccepted: " << m_trkCouldBeAccepted<< std::endl;
215
216 std::cout << "Dead Pixels: " << m_numPixelDeadSensor<< std::endl;
217 std::cout << "Dead SCT: " << m_numSCTDeadSensor << std::endl;
218 std::cout << "Pixel hit " << m_numPixelHits << std::endl;
219 std::cout << "Pixel holes " << m_numPixelHoles << std::endl;
220 std::cout << "SCT hits " << m_numSCTHits << std::endl;
221 std::cout << "SCT holes " << m_numSCTHoles << std::endl;
222
223 std::cout << "Unique hits " << m_numUnused << std::endl;
224 std::cout << "Unique SCT " << m_numSCT_Unused << std::endl;
225 std::cout << "Unique TRT " << m_numTRT_Unused << std::endl;
226 std::cout << "Unique Pix " << totalUniquePixelHits() << std::endl;
227 std::cout << "Pseudo " << m_numPseudo << std::endl;
228
229 std::cout << "SplitSharedPix " << m_numSplitSharedPix << std::endl;
230 std::cout << "SplitSharedSCT " << m_numSplitSharedSCT << std::endl;
231 std::cout << "Shared " << m_numShared << std::endl;
232 std::cout << "Shared SCT " << m_numSCT_Shared << std::endl;
233 std::cout << "Shared Pix " << totalSharedPixelHits() << std::endl;
234 std::cout << "Weighted shared " << m_numWeightedShared << std::endl;
235 };
236 };
237
239 bool operator() (const Trk::Track* x, const Trk::Track* y) const
240 {
241 if (ATH_UNLIKELY(!x and !y)) return false;
242 if (ATH_UNLIKELY(!x)) return true;
243 if (ATH_UNLIKELY(!y)) return false;
244 if (ATH_UNLIKELY(!x->trackParameters() and !y->trackParameters())) return false;
245 if (ATH_UNLIKELY(!x->trackParameters() || x->trackParameters()->size() <= 0) ) return true;
246 if (ATH_UNLIKELY(!y->trackParameters() || y->trackParameters()->size() <= 0) ) return false;
247 return std::fabs( (*x->trackParameters())[0]->parameters()[Trk::qOverP]) < std::fabs( (*y->trackParameters())[0]->parameters()[Trk::qOverP]) ;
248 }
249 };
250
252 {
253 unsigned int m_nTSoS{};
254 std::vector<int> m_type; // The type of TSOS
255 std::vector<int> m_detType; // The Detector type 1== Pixel, 11 = b-layer, 2 SCT, 3 TRT
256 std::vector<int> m_hitIsShared; // Number of tracks the hit is shared with
257 std::vector<float> m_splitProb1; // The Probablilty that the cluster on that surface is from 2 tracks
258 std::vector<float> m_splitProb2; // The Probablilty that the cluster on that surface is from 3 or more tracks
259 std::vector<const Trk::RIO_OnTrack*> m_RIO; // The cluster on track
260 std::multimap<const Trk::Track*, int, lessTrkTrack > m_overlappingTracks; // The tracks that overlap with the current track
261 std::multimap< int, const Trk::Track*> m_tracksSharingHit; // The tracks that overlap with the current track
262
264 {
265 std::cout << "WARNING DON'T USE THE DEFAULT CONSTRUCTOR OF tsosDetails" << std::endl;
266 };
267
268 TSoS_Details(unsigned int temp_nTSoS)
269 {
270 m_nTSoS = temp_nTSoS;
271 m_type.resize(m_nTSoS,OtherTsos);
272 m_detType.resize(m_nTSoS,-1.);
273 m_hitIsShared.resize(m_nTSoS, 0 );
274 m_splitProb1.resize(m_nTSoS,-1.f) ;
275 m_splitProb2.resize(m_nTSoS,-1.f) ;
276 m_RIO.resize(m_nTSoS,0);
277 };
278
279 int findIndexOfPreviousMeasurement( int currentIndex ) const
280 {
281 int indexPreviousMeasurement = currentIndex-1;
282 while(indexPreviousMeasurement >= 0){
283 if ( m_type[indexPreviousMeasurement] != OtherTsos ){
284 break;
285 } else {
286 --indexPreviousMeasurement;
287 }
288 } // end while
289 return indexPreviousMeasurement;
290 };
291 };
292
293 mutable std::mutex m_mutex;
294 struct CacheEntry {
295 EventContext::ContextEvt_t m_evt{EventContext::INVALID_CONTEXT_EVT};
296
297 // Max shared modules -- calulated from m_maxSharedModules
298 // 1 pixel hit is 1 module
299 // a double sided SCT hit (2 SCT hits) is 1 module
300 // so count by 2s for shared pixel hits and 1 per SCT (single sided hit) hit
302 // Min number of unique hits that are not already used on any other track
303 // but split hits can be used on multiple tracks and be considered unique
304 // - can change in ROI
306 // Min number of hits before we allow split sharing of hits -- can change if we are in ROI
308
309 std::vector<double> m_hadF;
310 std::vector<double> m_hadE;
311 std::vector<double> m_hadR;
312 std::vector<double> m_hadZ;
313
314 std::vector<double> m_emF;
315 std::vector<double> m_emE;
316 std::vector<double> m_emR;
317 std::vector<double> m_emZ;
318 };
319 mutable SG::SlotSpecificObj<CacheEntry> m_cache ATLAS_THREAD_SAFE; // Guarded by m_mutex
320
322 Trk::Track* createSubTrack( const std::vector<const Trk::TrackStateOnSurface*>& tsos, const Trk::Track* track ) const ;
323
325 void fillTrackDetails(const Trk::Track* ptrTrack,
326 Trk::ClusterSplitProbabilityContainer &splitProbContainer,
327 const Trk::PRDtoTrackMap &prd_to_track_map,
328 TrackHitDetails& trackHitDetails,
329 TSoS_Details& tsosDetails ) const;
330
338 const Trk::TrackScore score,
339 Trk::ClusterSplitProbabilityContainer &splitProbContainer,
340 Trk::PRDtoTrackMap &prd_to_track_map,
341 TrackHitDetails& trackHitDetails,
342 TSoS_Details& tsosDetails,
343 CacheEntry* ent,
344 int trackId) const;
345
350 bool performConversionCheck(const Trk::Track* ptrTrack,
351 Trk::PRDtoTrackMap &prd_to_track_map,
352 TrackHitDetails& trackHitDetails,
353 TSoS_Details& tsosDetails) const;
354
360 bool performHadDecayCheck(const Trk::Track* ptrTrack,
361 Trk::PRDtoTrackMap &prd_to_track_map,
362 TrackHitDetails& trackHitDetails,
363 TSoS_Details& tsosDetails) const;
364
369 void updateSharedForCollimated(TrackHitDetails& trackHitDetails,
370 TSoS_Details& tsosDetails) const;
371
372
375 Trk::ClusterSplitProbabilityContainer &clusterSplitProbMap) const;
376
378 bool inHadronicROI(const Trk::Track* ptrTrack) const;
379
382
384 bool isEmCaloCompatible(const Trk::TrackParameters& Tp) const;
385
387 static void newEvent(CacheEntry* ent) ;
388
393 bool checkOtherTracksValidity(TSoS_Details& tsosDetails,
394 int index,
395 Trk::ClusterSplitProbabilityContainer &splitProbContainer,
396 Trk::PRDtoTrackMap &prd_to_track_map,
397 int& maxiShared,
398 int& maxOtherNPixel,
399 bool& maxOtherHasIBL,
400 CacheEntry* ent) const;
401
403 std::pair< const Trk::TrackParameters* , const Trk::TrackParameters* >
404 getOverlapTrackParameters(int n, const Trk::Track* track1,
405 const Trk::Track* track2,
406 const Trk::PRDtoTrackMap &prd_to_track_map,
407 int splitSharedPix ) const;
408
412
414 bool clusCanBeSplit(float splitProb1, float splitProb2) const;
415 bool isTwoPartClus(float splitProb1, float splitProb2) const;
416 bool isMultiPartClus(float splitProb2) const;
417
418 static void rejectHitOverUse(TrackHitDetails& trackHitDetails, TSoS_Details& tsosDetails, int index) ;
419 static void rejectHit (TrackHitDetails& trackHitDetails, TSoS_Details& tsosDetails, int index) ;
420 static void rejectSharedHit (TrackHitDetails& trackHitDetails, TSoS_Details& tsosDetails, int index) ;
421 static void rejectSharedHitInvalid (TrackHitDetails& trackHitDetails, TSoS_Details& tsosDetails, int index) ;
422 static void sharedToSplitPix(TrackHitDetails& trackHitDetails, TSoS_Details& tsosDetails, int index) ;
423 static void addSharedHit (TrackHitDetails& trackHitDetails, TSoS_Details& tsosDetails, int index) ;
424 static void increaseSharedHitCounters(TrackHitDetails& trackHitDetails, bool isPix, bool isSCT) ;
425 static void decreaseSharedHitCounters(TrackHitDetails& trackHitDetails, bool isPix, bool isSCT) ;
426
428 PublicToolHandle<ITrtDriftCircleCutTool> m_selectortool{this, "DriftCircleCutTool", "InDet::InDetTrtDriftCircleCutTool"};
429
432
434 const SiliconID* m_detID{nullptr};
435
436 ToolHandle<Trk::IPRDtoTrackMapTool> m_assoTool
437 {this, "AssociationTool", "InDet::InDetPRDtoTrackMapToolGangedPixels" };
438
440 PublicToolHandle<Trk::ITrkObserverTool> m_observerTool{this, "ObserverTool", "", "track observer within ambiguity solver"};
441
443 IntegerProperty m_minHits{this, "minHits", 5, "Min Number of hits on track"};
444 IntegerProperty m_minTRT_Hits{this, "minTRTHits", 0, "Min Number of TRT hits on track"};
445 IntegerProperty m_maxSharedModules{this, "maxShared", 1, "Max number of shared modules"};
446 IntegerProperty m_maxTracksPerPRD{this, "maxTracksPerSharedPRD", 2, "Max number of tracks per hit. When NN is used, other flags set the limits."};
447 IntegerProperty m_minNotSharedHits{this, "minNotShared", 6, "Min number of non shared hits"};
448 FloatProperty m_minScoreShareTracks{this, "minScoreShareTracks", 0.0, "Min track score to alow it to share hits"};
449 BooleanProperty m_cosmics{this, "Cosmics", false, "Trying to reco cosmics?"};
450 BooleanProperty m_parameterization{this, "UseParameterization", true, "Use table of min number DCs"};
451 BooleanProperty m_doPixelClusterSplitting{this, "doPixelSplitting", false, "Split pixel clusters"};
452 FloatProperty m_sharedProbCut{this, "sharedProbCut", 0.3, "Min split prob to break a cluster into two parts"};
453 FloatProperty m_sharedProbCut2{this, "sharedProbCut2", 0.3, "Min split prob to break a clsuter into three parts"};
454
455 FloatProperty m_minTrackChi2ForSharedHits{this, "minTrackChi2ForSharedHits", 3, "Min track chi2 to split share hits"};
456 IntegerProperty m_minUniqueSCTHits{this, "minUniqueSCTHits", 2, "Min number of hits in the SCT that we need before we allow hit sharing in the SCT"};
457 IntegerProperty m_minSiHitsToAllowSplitting{this, "minSiHitsToAllowSplitting", 9, "Min number of hits before we allow split sharing of hits"};
458 IntegerProperty m_maxPixOnePartCluster{this, "maxPixOnePartCluster", 2, "Max number of tracks that can be associated to a 1 particle cluster"};
459 IntegerProperty m_maxPixTwoPartCluster{this, "maxPixTwoPartCluster", 2, "Max number of tracks that can be associated to a 2 particle cluster"};
460 IntegerProperty m_maxPixMultiCluster{this, "maxPixMultiCluster", 4, "Max number of tracks that can be associated to a >= 3 particle cluster"};
461 BooleanProperty m_shareSplitHits{this, "shareSplitHits", false, "Allow shared hits to be shared on 1 more track"};
462 IntegerProperty m_minPixHitAccepted{this, "minPixHitAccepted", 2, "Min number of pixel hits needed to be allowed to push accepted tracks over shared module limits"};
463
464 // ROI stuff
465 BooleanProperty m_useHClusSeed{this, "doHadCaloSeed", false};
466 BooleanProperty m_skipAmbiInROI{this, "doSkipAmbiInROI", false};
467 FloatProperty m_minPtSplit{this, "minPtSplit", 0.};
468 FloatProperty m_minPtBjetROI{this, "minPtBjetROI", 15000., "in MeV"};
469 FloatProperty m_phiWidth{this, "phiWidth", 0.2};
470 FloatProperty m_etaWidth{this, "etaWidth", 0.2};
472
473 BooleanProperty m_useEmClusSeed{this, "doEmCaloSeed", false};
474 FloatProperty m_phiWidthEm{this, "phiWidthEM", 0.05};
475 FloatProperty m_etaWidthEm{this, "etaWidthEM", 0.05};
476
478
479 //Track Pair Selection
480 BooleanProperty m_doPairSelection{this, "doPairSelection", true};
481 FloatProperty m_minPairTrackPt{this, "minPairTrackPt", 1000., "In MeV"};
482
483 BooleanProperty m_monitorTracks{this, "MonitorAmbiguitySolving", false, "to track observeration/monitoring (default is false)"};
484
485 };
486} // end of namespace
487
488
489#endif
#define ATH_UNLIKELY(x)
Maintain a set of objects, one per slot.
static Double_t Tp(Double_t *t, Double_t *par)
Property holding a SG store/key/clid from which a ReadHandle is made.
#define y
#define x
bool isEmCaloCompatible(const Trk::TrackParameters &Tp) const
Check if the cluster is compatible with a EM cluster.
static void rejectSharedHitInvalid(TrackHitDetails &trackHitDetails, TSoS_Details &tsosDetails, int index)
bool checkOtherTracksValidity(TSoS_Details &tsosDetails, int index, Trk::ClusterSplitProbabilityContainer &splitProbContainer, Trk::PRDtoTrackMap &prd_to_track_map, int &maxiShared, int &maxOtherNPixel, bool &maxOtherHasIBL, CacheEntry *ent) const
Returns true if accepted tracks remain about thresholds, false otherwise maxiShared = max number of s...
InDetDenseEnvAmbiTrackSelectionTool(const std::string &, const std::string &, const IInterface *)
void setPixelClusterSplitInformation(TSoS_Details &tsosDetails, Trk::ClusterSplitProbabilityContainer &clusterSplitProbMap) const
Update the pixel clusters split information.
SG::ReadHandleKey< ROIPhiRZContainer > m_inputHadClusterContainerName
virtual std::tuple< Trk::Track *, bool > getCleanedOutTrack(const Trk::Track *track, const Trk::TrackScore score, Trk::ClusterSplitProbabilityContainer &splitProbContainer, Trk::PRDtoTrackMap &prd_to_track_map, int trackId, int subtrackId) const override
Decide what to do with a candidate track.
bool isHadCaloCompatible(const Trk::TrackParameters &Tp) const
Check if the cluster is compatible with a hadronic cluster.
static void rejectHit(TrackHitDetails &trackHitDetails, TSoS_Details &tsosDetails, int index)
static void rejectSharedHit(TrackHitDetails &trackHitDetails, TSoS_Details &tsosDetails, int index)
SG::SlotSpecificObj< CacheEntry > m_cache ATLAS_THREAD_SAFE
bool inHadronicROI(const Trk::Track *ptrTrack) const
Does track pass criteria for hadronic ROI?
bool clusCanBeSplit(float splitProb1, float splitProb2) const
Simple helper functions to tell is cluster is split.
void updateSharedForCollimated(TrackHitDetails &trackHitDetails, TSoS_Details &tsosDetails) const
Handle update of the shared hit counts if either a conversion or a dense hadronic decay was identifie...
SG::ReadHandleKey< ROIPhiRZContainer > m_inputEmClusterContainerName
bool isTwoPartClus(float splitProb1, float splitProb2) const
static void sharedToSplitPix(TrackHitDetails &trackHitDetails, TSoS_Details &tsosDetails, int index)
ServiceHandle< IInDetEtaDependentCutsSvc > m_etaDependentCutsSvc
ITk eta-dependet cuts.
static void newEvent(CacheEntry *ent)
Fill hadronic & EM cluster map.
PublicToolHandle< Trk::ITrkObserverTool > m_observerTool
Observer tool This tool is used to observe the tracks and their 'score'.
PublicToolHandle< ITrtDriftCircleCutTool > m_selectortool
TRT minimum number of drift circles tool- returns allowed minimum number of TRT drift circles.
bool performConversionCheck(const Trk::Track *ptrTrack, Trk::PRDtoTrackMap &prd_to_track_map, TrackHitDetails &trackHitDetails, TSoS_Details &tsosDetails) const
Specific logic for identifing conversions with the goal of passing those tracks through to the final ...
std::pair< const Trk::TrackParameters *, const Trk::TrackParameters * > getOverlapTrackParameters(int n, const Trk::Track *track1, const Trk::Track *track2, const Trk::PRDtoTrackMap &prd_to_track_map, int splitSharedPix) const
Returns the Trackparameters of the two tracks on the n'th TrackStateOnSurface of the first track.
virtual ~InDetDenseEnvAmbiTrackSelectionTool()=default
default destructor
bool performHadDecayCheck(const Trk::Track *ptrTrack, Trk::PRDtoTrackMap &prd_to_track_map, TrackHitDetails &trackHitDetails, TSoS_Details &tsosDetails) const
Specific logic for identifing boosted light particle decays in jet topologies (tau and b),...
static void addSharedHit(TrackHitDetails &trackHitDetails, TSoS_Details &tsosDetails, int index)
void fillTrackDetails(const Trk::Track *ptrTrack, Trk::ClusterSplitProbabilityContainer &splitProbContainer, const Trk::PRDtoTrackMap &prd_to_track_map, TrackHitDetails &trackHitDetails, TSoS_Details &tsosDetails) const
Fill the two structs TrackHitDetails & TSoS_Details full of information.
static void increaseSharedHitCounters(TrackHitDetails &trackHitDetails, bool isPix, bool isSCT)
Trk::Track * createSubTrack(const std::vector< const Trk::TrackStateOnSurface * > &tsos, const Trk::Track *track) const
method to create a new track from a vector of TSOS's
virtual StatusCode initialize() override
standard Athena-Algorithm method
static void rejectHitOverUse(TrackHitDetails &trackHitDetails, TSoS_Details &tsosDetails, int index)
static void decreaseSharedHitCounters(TrackHitDetails &trackHitDetails, bool isPix, bool isSCT)
bool isNearbyTrackCandidate(const Trk::TrackParameters *paraA, const Trk::TrackParameters *paraB) const
Check if two sets of track paremeters are compatible with being from a the same low mass particle dec...
void decideWhichHitsToKeep(const Trk::Track *, const Trk::TrackScore score, Trk::ClusterSplitProbabilityContainer &splitProbContainer, Trk::PRDtoTrackMap &prd_to_track_map, TrackHitDetails &trackHitDetails, TSoS_Details &tsosDetails, CacheEntry *ent, int trackId) const
Determine which hits to keep on this track Look at the hits on track and decided if they should be ke...
Property holding a SG store/key/clid from which a ReadHandle is made.
Maintain a set of objects, one per slot.
This is an Identifier helper class for both the Pixel and SCT subdetectors.
Definition SiliconID.h:42
Container to associate Cluster with cluster splitting probabilities.
Class to handle RIO On Tracks ROT) for InDet and Muons, it inherits from the common MeasurementBase.
Definition RIO_OnTrack.h:70
represents the track state (measurement, material, fit parameters and quality) at a surface.
Primary Vertex Finder.
Ensure that the ATLAS eigen extensions are properly loaded.
float TrackScore
Definition TrackScore.h:10
@ qOverP
perigee
Definition ParamDefs.h:67
ParametersBase< TrackParametersDim, Charged > TrackParameters
Definition index.py:1
std::multimap< const Trk::Track *, int, lessTrkTrack > m_overlappingTracks
bool operator()(const Trk::Track *x, const Trk::Track *y) const