ATLAS Offline Software
Loading...
Searching...
No Matches
TrackStatHelper.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
6// StatBox.h
7// Authors: Sven Vahsen
9
10// Private utility class used by IDTrackStat to do track counting
11// Since this StatBox is likely to change to without notification
12// as needed, I would not recommend anyone else to use it
13
14#ifndef INDETRECSTATISTICS_TrackStatHelper_H
15#define INDETRECSTATISTICS_TrackStatHelper_H
16
17#include "TrkTrack/TrackInfo.h"
18#include "TrkTrackSummary/TrackSummary.h" // needed by Trk::numberOfTrackSummaryTypes
24
25#include <vector>
26#include <string>
27#include <map>
28#include <atomic>
29
30// forward declarations:
31class PixelID;
32class SCT_ID;
33class AtlasDetectorID;
34class Track;
35
36namespace InDet {
37
40 TRACK_SECONDARY, //LT 06.21
73
82
83 struct cuts {
100 float minPt;
102 maxEtaBarrel =0;
104 maxEtaEndcap =0;
105 minEtaFORWARD = 2.5;
106 maxEtaFORWARD=4.2;
107 fakeTrackCut=0;
108 fakeTrackCut2 =0;
109 matchTrackCut =0;
110 maxRStartPrimary = 25.0;
111 maxRStartSecondary= 360.0;
112 maxZStartPrimary = 200.0;
113 maxZStartSecondary=2000.0;
114 minREndPrimary = 400.0;
115 minREndSecondary =1000.0;
116 minZEndPrimary =2300.0;
117 minZEndSecondary =3200.0;
118 minPt = 1000;
119 }
120
121
122
123 };
124
125
132
133 public:
136 TrackStatHelper (const std::string&, const std::string&, bool careAboutTruth = true);
138 void SetCuts(const struct cuts&);
140 void addEvent (const EventContext& ctx,
141 const TrackCollection *,
142 std::vector<const Trk::Track *> &,
143 const std::vector <std::pair<HepMC::ConstGenParticlePtr,int> > &,
144 const TrackTruthCollection *,
145 const AtlasDetectorID * const,
146 const PixelID *,
147 const SCT_ID *,
149 bool,
150 const unsigned int *,
151 const unsigned int *) const;
153 void reset ();
155 void print (MsgStream &out) const;
157 void printRegion1(MsgStream &out, enum eta_region) const;
158 void printRegion2(MsgStream &out, enum eta_region, float denominator) const;
160 void printSecondary (MsgStream &out) const;
162 void printRegionSecondary(MsgStream &out, enum eta_region, float denominator) const;
163
164
166 bool printTrackSummaryRegion (MsgStream &out, enum track_types, enum eta_region) const;
168 void printTrackSummaryAverage(MsgStream &out, enum track_types, enum eta_region, int summary_type) const;
170 bool PassTrackCuts(const Trk::TrackParameters *para) const;
172 int ClassifyParticle( const HepMC::ConstGenParticlePtr& particle, const double prob) const;
173
174
175 static std::string getSummaryTypeHeader();
176
178 const std::string &key() const { return m_TrackCollectionKey; }
180 const std::string &Truthkey() const { return m_TrackTruthCollectionKey; }
181
182 private:
183
186
188 mutable std::atomic<long> m_events {} ;
189
190 template <int N_Categories, int N_Types, int N_Regions, typename T_Int=long>
191 struct Counter {
192 Counter() { reset(); }
193 void reset() {
194 for (unsigned int cat_i=0; cat_i < N_Categories; ++cat_i ) {
195 for (unsigned int type_i=0; type_i < N_Types; ++type_i) {
196 for (unsigned int eta_i=0; eta_i < N_Regions; ++eta_i) {
197 m_counter[cat_i][type_i][eta_i]=0;
198 }
199 }
200 }
201 }
202 template <typename T_IntB>
204 for (unsigned int cat_i=0; cat_i < N_Categories; ++cat_i ) {
205 for (unsigned int type_i=0; type_i < N_Types; ++type_i) {
206 for (unsigned int eta_i=0; eta_i < N_Regions; ++eta_i) {
207 m_counter[cat_i][type_i][eta_i] += a.m_counter[cat_i][type_i][eta_i];
208 }
209 }
210 }
211 return *this;
212 }
213 T_Int m_counter[N_Categories][N_Types][N_Regions];
214 };
215
216 template <int N_Categories, int N_Types, int N_Regions, int N_SubCategories, typename T_Int=long>
217 struct Counter4D {
219 void reset() {
220 for (unsigned int cat_i=0; cat_i < N_Categories; ++cat_i ) {
221 for (unsigned int type_i=0; type_i < N_Types; ++type_i) {
222 for (unsigned int eta_i=0; eta_i < N_Regions; ++eta_i) {
223 for (unsigned int sub_i=0; sub_i < N_SubCategories; ++sub_i) {
224 m_counter[cat_i][type_i][eta_i][sub_i]=0;
225 }
226 }
227 }
228 }
229 }
230 template <typename T_IntB>
232 for (unsigned int cat_i=0; cat_i < N_Categories; ++cat_i ) {
233 for (unsigned int type_i=0; type_i < N_Types; ++type_i) {
234 for (unsigned int eta_i=0; eta_i < N_Regions; ++eta_i) {
235 for (unsigned int sub_i=0; sub_i < N_SubCategories; ++sub_i) {
236 m_counter[cat_i][type_i][eta_i][sub_i] += a.m_counter[cat_i][type_i][eta_i][sub_i];
237 }
238 }
239 }
240 }
241 return *this;
242 }
243 T_Int m_counter[N_Categories][N_Types][N_Regions][N_SubCategories];
244 };
245
250
254
270
271 // The ETrackSummaryTypes should be synchronised with the arrays summaryTypes and
272 // summaryTypeName i.e. matching order and number of elements. Of these enums only
273 // two are used: kNSummaryTypes which defines the summary statistics array sizes;
274 // kNumberOfPixelHits, which is used to get counts for some denominators.
292 static const char *const s_summaryTypeName [kNSummaryTypes];
294
296 N_TRACKTYPES,
297 N_ETAREGIONS,
299 int>;
301 N_TRACKTYPES,
302 N_ETAREGIONS,
304 std::atomic<long> >;
306
307
308 void setSummaryStat(track_types track_i,
309 eta_region region_i,
310 const Trk::TrackSummary *summary,
311 TrackSummaryCounter &trackSummarySum) const
312 {
313 if (summary) {
314 for (int stype=0; stype < kNSummaryTypes; stype++) {
315 int value = summary->get(s_summaryTypes[stype]);
316 //value is -1 if undefined
317 if (value>0) {
318 trackSummarySum.m_counter[kTrackSummarySum][track_i][region_i][stype] += value;
319 trackSummarySum.m_counter[kNTrackSummaryOK][track_i][region_i][stype] ++;
320 }
321 else {
322 trackSummarySum.m_counter[kNTrackSummaryBAD][track_i][region_i][stype] ++;
323 }
324 }
325 }
326 }
327
328 mutable std::atomic<bool> m_truthMissing;
330 struct cuts m_cuts;
331
332 mutable std::mutex m_authorMutex;
333 mutable std::bitset<Trk::TrackInfo::NumberOfTrackRecoInfo> m_recoInfo ATLAS_THREAD_SAFE;
334 mutable std::bitset<Trk::TrackInfo::NumberOfTrackProperties> m_patternProperties ATLAS_THREAD_SAFE;
335
336 typedef std::multimap<HepMcParticleLink,float> recoToTruthMap;
337 };
338
339
340
341} // close of namespace
342
343#endif // INDETRECSTATISTICS_TrackStatHelper_H
static Double_t a
DataVector< Trk::Track > TrackCollection
This typedef represents a collection of Trk::Track objects.
This class provides an interface to generate or decode an identifier for the upper levels of the dete...
bool printTrackSummaryRegion(MsgStream &out, enum track_types, enum eta_region) const
Sets up detailed statistics part of table, calls printTrackSummaryRegion.
bool PassTrackCuts(const Trk::TrackParameters *para) const
defines 'good' reco tracks
void printRegionSecondary(MsgStream &out, enum eta_region, float denominator) const
Prints ntracks per event,efficiencies,fake rates, and general hit information for given eta region.
void printRegion1(MsgStream &out, enum eta_region) const
Prints ntracks per event,efficiencies,fake rates, and general hit information for given eta region.
std::atomic< bool > m_author_found[Trk::TrackInfo::NumberOfTrackFitters]
Number of tracking authors found.
Counter4D< kNTrackSummaryCounter, N_TRACKTYPES, N_ETAREGIONS, kNSummaryTypes, int > TrackSummaryCounter
type statistics.
TrackStatHelper(const std::string &, const std::string &, bool careAboutTruth=true)
Constructor.
void printSecondary(MsgStream &out) const
Prints all of the statistics information, calls printRegion, printTrackSummaryRegion,...
std::string m_TrackTruthCollectionKey
StoreGate Track Truth Collection Key.
std::string m_TrackCollectionKey
StoreGate Track Collection Key.
TracksCounterAtomic m_tracks ATLAS_THREAD_SAFE
void printTrackSummaryAverage(MsgStream &out, enum track_types, enum eta_region, int summary_type) const
Prints information from TrackSummaryTool.
static std::string getSummaryTypeHeader()
static const char *const s_summaryTypeName[kNSummaryTypes]
table column labels for summary
std::multimap< HepMcParticleLink, float > recoToTruthMap
map containing reco track and matched truth track barcode
const std::string & Truthkey() const
Returns Truth TrackCollection Key.
Counter< kNTracksCounter, N_TRACKTYPES, N_ETAREGIONS, int > TracksCounter
void print(MsgStream &out) const
Prints all of the statistics information, calls printRegion, printTrackSummaryRegion,...
@ kTracks_rec
number of reconstructed tracks for a given type and eta region
@ kTracks_gen
number of generated tracks for a given type and eta region, looping over genevents to include possibl...
@ kTracks_gen_signal
number of generated tracks for a given type and eta region, just from first genevent
Counter< kNTracksCounter, N_TRACKTYPES, N_ETAREGIONS, std::atomic< long > > TracksCounterAtomic
static const Trk::SummaryType s_summaryTypes[kNSummaryTypes]
summary types for which statistics are gathered
void reset()
Resets the track collection information, called in the constructor.
void setSummaryStat(track_types track_i, eta_region region_i, const Trk::TrackSummary *summary, TrackSummaryCounter &trackSummarySum) const
void addEvent(const EventContext &ctx, const TrackCollection *, std::vector< const Trk::Track * > &, const std::vector< std::pair< HepMC::ConstGenParticlePtr, int > > &, const TrackTruthCollection *, const AtlasDetectorID *const, const PixelID *, const SCT_ID *, const Trk::IExtendedTrackSummaryTool *, bool, const unsigned int *, const unsigned int *) const
Adds hit, track and matching information for each event.
@ kNTrackSummaryBAD
Number of tracks with track summary bad for given type,eta,summary type.
@ kNTrackSummaryOK
Number of tracks with track summary OK for given type,eta,summary type.
@ kTrackSummarySum
Track Summary Values for each track type, region and summary type.
const std::string & key() const
Returns TrackCollection Key.
Counter< kNHitsCounter, N_HITTYPES, N_ETAREGIONS, std::atomic< long > > HitsCounterAtomic
Counter4D< kNTrackSummaryCounter, N_TRACKTYPES, N_ETAREGIONS, kNSummaryTypes, std::atomic< long > > TrackSummaryCounterAtomic
@ kHits_sec
number of hits from secondary tracks for a given type and eta region
@ kHits_rec
number of reconstructed hits for a given type and eta region
@ kHits_pri
number of hits from primary tracks for a given type and eta region
int ClassifyParticle(const HepMC::ConstGenParticlePtr &particle, const double prob) const
classifies gen particle as primary, secondary or truncated
Counter< kNHitsCounter, N_HITTYPES, N_ETAREGIONS, int > HitsCounter
std::atomic< long > m_events
Number of events.
void printRegion2(MsgStream &out, enum eta_region, float denominator) const
std::atomic< bool > m_truthMissing
Flag for if track truth is missing.
void SetCuts(const struct cuts &)
Sets the cuts such as the eta regions (barrel, transition,endcap) and the hit fraction fake cuts and ...
@ TRACK_SECONDARY
@ TRACK_NOHEPMCPARTICLELINK
@ TRACK_LOWTRUTHPROB2_SIGNAL
@ TRACK_LOWTRUTHPROB2
@ TRACK_TRUNCATED
@ TRACK_MATCHED_PRIMARY
@ TRACK_LOWTRUTHPROB_SIGNAL
@ TRACK_MULTMATCH_SECONDARY
@ TRACK_MULTMATCH_PRIMARY
@ TRACK_MULTMATCH
@ TRACK_MATCHED_SIGNAL
@ TRACK_MATCHED_SECONDARY
@ TRACK_LOWTRUTHPROB
This is an Identifier helper class for the Pixel subdetector.
Definition PixelID.h:69
This is an Identifier helper class for the SCT subdetector.
Definition SCT_ID.h:68
Interface for condensing Trk::Track properties and associated hits to a (non-fittable) foot print,...
A summary of the information contained by a track.
const GenParticle * ConstGenParticlePtr
Definition GenParticle.h:38
Primary Vertex Finder.
ParametersBase< TrackParametersDim, Charged > TrackParameters
SummaryType
enumerates the different types of information stored in Summary.
Counter4D & operator+=(const Counter4D< N_Categories, N_Types, N_Regions, N_SubCategories, T_IntB > &a)
Counter & operator+=(const Counter< N_Categories, N_Types, N_Regions, T_IntB > &a)
float maxRStartSecondary
float minZEndSecondary
float fakeTrackCut
fraction of hits per track that come from single matched truth track.
float maxRStartPrimary
float matchTrackCut
Truth probability has to be greater than this for track to be considered matched.
float minREndSecondary
float maxEtaTransition
Maxiumu eta for transition region.
float maxZStartPrimary
float fakeTrackCut2
2nd value for fraction of hits per track that come from single matched truth track.
float maxEtaEndcap
Maximum eta for endcap.
float maxZStartSecondary