ATLAS Offline Software
Loading...
Searching...
No Matches
TrackStatHelper.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
6// TrackStatHelper.cxx
7// Authors: Sven Vahsen
9
10// Private utility class used by IDTrackStat to do track and hit counting
11// Since I'm likely to changes the interface as needed,
12// I would not recommend anyone to use this class outside of IDTrackStat
13//
14// to-do-list
15// o remove TrackCollection argument to addEvent. Is this really needed for the truthmatching?
16// o add DR cut to the fake track cut
17//
18
22#include "AtlasHepMC/GenEvent.h"
23#include "CLHEP/Geometry/Point3D.h"
30#include "TrkTrack/Track.h"
33#include <cstring>
34#include <iomanip>
35#include <sstream>
36#include <utility>
37#include <vector>
38
39namespace Trk {
40
42 public:
43 static void extract(std::vector<const RIO_OnTrack*>& rots, const std::vector<const MeasurementBase*>& measurements);
44 };
45
46}
47
48
49// summary types for which statistics are gathered.
50// should be synchronised with ETrackSummaryTypes
64};
65
66// Table column labels. should be synchronised with ETrackSummaryTypes
67const char * const InDet::TrackStatHelper::s_summaryTypeName[kNSummaryTypes] = {
68 "blay",
69 "pix",
70 "hole",
71 "gang",
72 "SCT",
73 "hole",
74 "DHole",
75 "TRT",
76 "outl",
77 "TRHi",
78 "outl",
79 "alloutl"};
80
82 std::stringstream out;
83 for (unsigned int stype_i=0; stype_i < kNSummaryTypes; ++stype_i ) {
84 out << std::setw(std::max(6,static_cast<int>(strlen(s_summaryTypeName[stype_i])+1))) << s_summaryTypeName[stype_i];
85 }
86 return out.str();
87}
88
89
90static const std::string track_types_string[InDet::N_TRACKTYPES]= {"all",
91 "primary",
92 "secondary",
93 "truncated",
94 "lowPb1",
95 "lowPb2",
96 "noLink",
97 "matched",
98 "matchedPrimary",
99 "matchedSecondary",
100 "multipleMatches"
101 };//LT added 06.21
102
103InDet::TrackStatHelper::TrackStatHelper(const std::string& TrackCollectionKey, const std::string& TrackTruthCollectionKey, bool careAboutTruth):
104 m_TrackCollectionKey (std::move(TrackCollectionKey)),
105 m_TrackTruthCollectionKey (std::move(TrackTruthCollectionKey)),
106 m_truthMissing (false),
107 m_careAboutTruth (careAboutTruth)
108{
109 this->reset();
110}
111
112
114{
115 m_cuts = ct;
116
117}
118
119void InDet::TrackStatHelper::addEvent(const EventContext& ctx,
120 const TrackCollection * recTracks,
121 std::vector <const Trk::Track *> & rec,
122 const std::vector <std::pair<HepMC::ConstGenParticlePtr,int> > & gen,
123 const TrackTruthCollection * truthMap,
124 const AtlasDetectorID * const idHelper,
125 const PixelID * pixelID,
126 const SCT_ID * sctID,
127 const Trk::IExtendedTrackSummaryTool* trkSummaryTool,
128 bool useTrackSummary,
129 const unsigned int * inTimeStart,
130 const unsigned int * inTimeEnd) const
131{
132
133 recoToTruthMap rttMap;
134
135 recoToTruthMap::const_iterator imap;
136
137 m_events ++;
138
139 float Eta = 0;
140 InDet::eta_region Region = ETA_ALL;
141 int Author = 0;
142 int recoClassification = 0;
143 bool truth_missing=false;
144
145 TracksCounter tracks;
146 HitsCounter hits;
147 TrackSummaryCounter trackSummarySum;
148 bool author_found [Trk::TrackInfo::NumberOfTrackFitters];
149 std::bitset<Trk::TrackInfo::NumberOfTrackRecoInfo> reco_info;
150 std::bitset<Trk::TrackInfo::NumberOfTrackProperties> pattern_properties;
151 for (unsigned int i=0; i < Trk::TrackInfo::NumberOfTrackFitters; ++i) { author_found[i]=false; }
152
153 // ------------ reconstructed tracks -----------------
154 for (const Trk::Track* track : rec) {
155
156 const Trk::TrackParameters* para = track->trackParameters()->front();
157 recoClassification = -999;
158 int thisEventIndex = -999;
159
160 Author = track->info().trackFitter();
161 if (Author > 0 && Author < Trk::TrackInfo::NumberOfTrackFitters){
162 author_found[Author] = true;
163 reco_info |= track->info().patternRecognition();
164 pattern_properties |= track->info().properties();
165 }
166 else {
167 //FIX author_problem++;
168 }
169 if (!para) {
170 Region = ETA_UNKNOWN;
171 }
172 else {
173 Eta = std::abs(para->eta());
174 if (Eta < m_cuts.maxEtaBarrel) Region = ETA_BARREL;
175 else if (Eta < m_cuts.maxEtaTransition) Region = ETA_TRANSITION;
176 else if (Eta < m_cuts.maxEtaEndcap) Region = ETA_ENDCAP;
177 else if ((Eta > m_cuts.minEtaFORWARD) && (Eta < m_cuts.maxEtaFORWARD)) Region = ETA_FORWARD;
178 else Region = ETA_OUTSIDE;
179 }
180
181 if ( !(Region==ETA_BARREL || Region==ETA_TRANSITION || Region==ETA_ENDCAP || Region==ETA_FORWARD) )
182 {
183 continue; // Only want to tally tracks that are within the detector volume.
184 }
185
186 TrackTruthCollection::const_iterator found;
187 TrackTruth trtruth;
188 double trprob = 0;
189
190 if (!truthMap) {
191 // no truthmap
192 truth_missing=true;
193 }
194 else {
196 tracklink.setElement(track);
197 tracklink.setStorableObject(*recTracks);
198 const ElementLink<TrackCollection> tracklink2=tracklink;
199
200 found = truthMap->find(tracklink2);
201 if(!(found == truthMap->end())){
202 trtruth=found->second;
203 HepMcParticleLink hmpl = trtruth.particleLink();
204 thisEventIndex = hmpl.eventIndex();
205 rttMap.insert(std::pair<HepMcParticleLink,float>(hmpl,trtruth.probability()));
206 // ME: remember prob
207 trprob = trtruth.probability();
208 }
209
210 }
211
212
213 if(!PassTrackCuts(para)) continue;
214
215 tracks.m_counter[kTracks_rec][TRACK_ALL][Region]++;
216 tracks.m_counter[kTracks_rec][TRACK_ALL][ETA_ALL]++;
217 // signal only tracks for the denominator in signal efficiencies
218 if (thisEventIndex==0)
219 {
220 tracks.m_counter[kTracks_rec][TRACK_ALL_SIGNAL][Region]++;
221 tracks.m_counter[kTracks_rec][TRACK_ALL_SIGNAL][ETA_ALL]++;
222 }
223
224 // process track summary
225
226 std::unique_ptr<Trk::TrackSummary> cleanup;
227 const Trk::TrackSummary* summary = track->trackSummary();
228
229 if (useTrackSummary) {
230 if (!track->trackSummary()) {
231 cleanup = trkSummaryTool->summary(ctx, *track);
232 summary=cleanup.get();
233 }
234
235 setSummaryStat(TRACK_ALL, Region, summary, trackSummarySum);
236 }
237
238 // check if current reco track has more than a specified fraction of hits
239 // from a generated track. [If not, it will be counted as a fake
240 // track in the statistics report.]
241
242 if (!truthMap) {
243 truth_missing=true;
244 }
245 else {// no truthmap
246
247 if (found == truthMap->end()) {
248 // m_truthMissing=true;
249 // no truth might happen with new truth !
250 tracks.m_counter[kTracks_rec][TRACK_NOHEPMCPARTICLELINK][Region]++;
251 tracks.m_counter[kTracks_rec][TRACK_NOHEPMCPARTICLELINK][ETA_ALL]++;
252 setSummaryStat(TRACK_NOHEPMCPARTICLELINK, Region, summary, trackSummarySum);
253 }
254 else{ // no link
255
256 // ME : change logic, secondaries from G4 processes are truncated
257 const HepMcParticleLink& HMPL=trtruth.particleLink();
258 if (! HMPL.isValid()) {
259 tracks.m_counter[kTracks_rec][TRACK_NOHEPMCPARTICLELINK][Region]++;
260 tracks.m_counter[kTracks_rec][TRACK_NOHEPMCPARTICLELINK][ETA_ALL]++;
261 setSummaryStat(TRACK_NOHEPMCPARTICLELINK, Region, summary, trackSummarySum);
262 }
263 else {
264 //classify track as coming from primary, secondary or truncated gen particle
265#ifdef HEPMC3
266 HepMC::ConstGenParticlePtr particle = HMPL.scptr();
267#else
268 const HepMC::GenParticle *particle = HMPL.cptr();
269#endif
270 recoClassification = ClassifyParticle(particle, trprob);
271
272 if (trprob < m_cuts.fakeTrackCut)
273 {
274 tracks.m_counter[kTracks_rec][TRACK_LOWTRUTHPROB][Region]++;
275 tracks.m_counter[kTracks_rec][TRACK_LOWTRUTHPROB][ETA_ALL]++;
276 if (thisEventIndex==0)
277 {
278 tracks.m_counter[kTracks_rec][TRACK_LOWTRUTHPROB_SIGNAL][Region]++;
279 tracks.m_counter[kTracks_rec][TRACK_LOWTRUTHPROB_SIGNAL][ETA_ALL]++;
280 }
281 setSummaryStat(TRACK_LOWTRUTHPROB, Region, summary, trackSummarySum);
282 }
283 if (trprob < m_cuts.fakeTrackCut2) {
284 tracks.m_counter[kTracks_rec][TRACK_LOWTRUTHPROB2][Region]++;
285 tracks.m_counter[kTracks_rec][TRACK_LOWTRUTHPROB2][ETA_ALL]++;
286 if (thisEventIndex==0)
287 {
288 tracks.m_counter[kTracks_rec][TRACK_LOWTRUTHPROB2_SIGNAL][Region]++;
289 tracks.m_counter[kTracks_rec][TRACK_LOWTRUTHPROB2_SIGNAL][ETA_ALL]++;
290 }
291 setSummaryStat(TRACK_LOWTRUTHPROB2, Region, summary, trackSummarySum);
292 }
293 } // end else on !isVadid
294 }
295 }
296
297 // ------------------ hits on reconstructed tracks ---------------
298
299 EHitsCounter part_type = ( recoClassification==TRACK_PRIMARY
300 ? kHits_pri
301 : ( recoClassification==TRACK_SECONDARY
302 ? kHits_sec
303 : kNHitsCounter));
304 for (const Trk::TrackStateOnSurface* hit : *track->trackStateOnSurfaces()) {
305 if(hit){
306 const Trk::MeasurementBase* mesh =hit->measurementOnTrack();
307 if (mesh) {
308 const Trk::RIO_OnTrack* rio{} ;
310 assert(dynamic_cast<const Trk::RIO_OnTrack*>(mesh)!=nullptr);
311 rio = static_cast<const Trk::RIO_OnTrack*>(mesh);
312 } else {
313 // try CompetingROT:
315 assert(dynamic_cast<const Trk::CompetingRIOsOnTrack*>(mesh) !=nullptr);
316 auto comprot = static_cast<const Trk::CompetingRIOsOnTrack*>(mesh);
317 rio = &comprot->rioOnTrack(comprot->indexOfMaxAssignProb());
318 }
319 }
320
321
322 if (rio) {
323 // skip outliers
325 continue;
326
327 hits.m_counter[kHits_rec][HIT_ALL][Region]++;
328 hits.m_counter[kHits_rec][HIT_ALL][ETA_ALL]++;
329 if(part_type<kNHitsCounter){
330 hits.m_counter[part_type][HIT_ALL][Region]++;
331 hits.m_counter[part_type][HIT_ALL][ETA_ALL]++;
332 }
333
334 Identifier id = rio->identify();
335 int HitDet = HIT_UNKNOWN;
336 int HitLayer = N_HITTYPES;
337 bool part_type_for_all=true;
338 if (idHelper->is_trt(id)){
339 HitDet=HIT_TRT_ALL;
340 }
341 else if (idHelper->is_sct(id)){
342 HitDet=HIT_SCT_ALL;
343 if (sctID) {
344 switch (sctID->layer_disk(id)) {
345 case 0: HitLayer = HIT_SCT1; break;
346 case 1: HitLayer = HIT_SCT2; break;
347 case 2: HitLayer = HIT_SCT3; break;
348 case 3: HitLayer = HIT_SCT4; break;
349 case 4: case 5: case 6: case 7: case 8: HitLayer = HIT_SCT5TO9; break;
350 default: HitLayer = HIT_UNKNOWN; break;
351 }
352 }
353 }
354 else if (idHelper->is_pixel(id)){
355 part_type_for_all=false;
356 HitDet = HIT_PIXEL_ALL;
357 if (pixelID) {
358 switch (pixelID->layer_disk(id)) {
359 case 0: HitLayer = HIT_PIX1; break;
360 case 1: HitLayer = HIT_PIX2; break;
361 case 2: HitLayer = HIT_PIX3; break;
362 default: HitLayer = HIT_UNKNOWN;
363 }
364 }
365 }
366
367 hits.m_counter[kHits_rec][HitDet][Region]++;
368 hits.m_counter[kHits_rec][HitDet][ETA_ALL]++;
369 if(part_type_for_all && part_type<kNHitsCounter){
370 hits.m_counter[part_type][HitDet][Region]++;
371 hits.m_counter[part_type][HitDet][ETA_ALL]++;
372 }
373
374 if (HitLayer<N_HITTYPES) {
375 hits.m_counter[kHits_rec][HitLayer][Region]++;
376 hits.m_counter[kHits_rec][HitLayer][ETA_ALL]++;
377 if(part_type<kNHitsCounter){
378 hits.m_counter[part_type][HitLayer][Region]++;
379 hits.m_counter[part_type][HitLayer][ETA_ALL]++;
380 }
381 }
382 }
383 }
384 }
385 }
386 }
387
388 // ------------------------- generated tracks, including pileup -----------------------------
389
390 Eta = 0;
391 Region = ETA_ALL;
392 int classification=-999;
393 for (auto truth = gen.begin(); truth != gen.end(); ++truth) {
394 classification=-999;
395 bool inTimePileup = truth->second == 0 || (truth->second >= (int)*inTimeStart && truth->second <= (int)*inTimeEnd);
396
397 auto particle = truth->first;
398
399 //determine eta region
400 Eta = std::abs(particle->momentum().pseudoRapidity());
401 if (Eta < m_cuts.maxEtaBarrel) Region = ETA_BARREL;
402 else if (Eta < m_cuts.maxEtaTransition) Region = ETA_TRANSITION;
403 else if (Eta < m_cuts.maxEtaEndcap) Region = ETA_ENDCAP;
404 else if ((Eta > m_cuts.minEtaFORWARD) && (Eta < m_cuts.maxEtaFORWARD)) Region = ETA_FORWARD;
405 else Region = ETA_OUTSIDE;
406 if ( !(Region==ETA_BARREL || Region==ETA_TRANSITION || Region==ETA_ENDCAP || Region==ETA_FORWARD) )
407 {
408 continue; // Only want to tally tracks that are within the detector volume.
409 }
410 tracks.m_counter[kTracks_gen][TRACK_ALL][Region]++;
411 tracks.m_counter[kTracks_gen][TRACK_ALL][ETA_ALL] ++;
412 if (inTimePileup) tracks.m_counter[kTracks_gen][TRACK_ALL_SIGNAL][ETA_ALL]++; // "SIGNAL" is misnomer here
413 if (inTimePileup) tracks.m_counter[kTracks_gen][TRACK_ALL_SIGNAL][Region]++;
414
415 //classify gen as primary, secondary or truncated
416 classification = ClassifyParticle(particle,1.);
417
418 if(classification==TRACK_PRIMARY){
419 tracks.m_counter[kTracks_gen][TRACK_PRIMARY][ETA_ALL] ++;
420 tracks.m_counter[kTracks_gen][TRACK_PRIMARY][Region]++;
421 }
422 if(classification==TRACK_SECONDARY){
423 tracks.m_counter[kTracks_gen][TRACK_SECONDARY][ETA_ALL] ++;
424 tracks.m_counter[kTracks_gen][TRACK_SECONDARY][Region]++;
425 }
426
427 //see if gen track has at least 1 matching reco track with high enough probability
428 bool matched = false;
429 int nmatched = 0;
430
431 HepMcParticleLink hmpl2(particle,particle->parent_event()->event_number(),HepMcParticleLink::IS_EVENTNUM);
432 recoToTruthMap::iterator rttIter=rttMap.find(hmpl2);
433 if(rttIter != rttMap.end()){
434 for(imap = rttMap.lower_bound(hmpl2); imap !=rttMap.upper_bound(hmpl2); ++imap){
435 if(imap->second > m_cuts.matchTrackCut){
436 matched = true;
437 nmatched++;
438 }
439 }
440 }
441 if (matched) {
442 tracks.m_counter[kTracks_gen][TRACK_MATCHED][Region]++;
443 tracks.m_counter[kTracks_gen][TRACK_MATCHED][ETA_ALL]++;
444 if (inTimePileup)
445 {
446 tracks.m_counter[kTracks_gen][TRACK_MATCHED_SIGNAL][Region]++;
447 tracks.m_counter[kTracks_gen][TRACK_MATCHED_SIGNAL][ETA_ALL]++;
448 }
449 if(nmatched > 1){
450 tracks.m_counter[kTracks_gen][TRACK_MULTMATCH][Region] += nmatched-1;
451 tracks.m_counter[kTracks_gen][TRACK_MULTMATCH][ETA_ALL]+= nmatched-1;
452 }
453 if(classification==TRACK_PRIMARY){
454 tracks.m_counter[kTracks_gen][TRACK_MATCHED_PRIMARY][Region]++;
455 tracks.m_counter[kTracks_gen][TRACK_MATCHED_PRIMARY][ETA_ALL]++;
456 if(nmatched > 1){
457 tracks.m_counter[kTracks_gen][TRACK_MULTMATCH_PRIMARY][Region] += nmatched-1;
458 tracks.m_counter[kTracks_gen][TRACK_MULTMATCH_PRIMARY][ETA_ALL] += nmatched-1;
459 }
460 }
461 if(classification==TRACK_SECONDARY){
462 tracks.m_counter[kTracks_gen][TRACK_MATCHED_SECONDARY][Region]++;
463 tracks.m_counter[kTracks_gen][TRACK_MATCHED_SECONDARY][ETA_ALL]++;
464 if(nmatched > 1){
465 tracks.m_counter[kTracks_gen][TRACK_MULTMATCH_SECONDARY][Region] += nmatched-1;
466 tracks.m_counter[kTracks_gen][TRACK_MULTMATCH_SECONDARY][ETA_ALL] += nmatched-1;
467 }
468 }
469 }
470 }
471
472 // ------------------------- generated tracks, just signal event -----------------------------
473 // For efficiencies.
474
475 Eta = 0;
476 Region = ETA_ALL;
477 classification=-999;
478
479 for (auto truth = gen.begin(); truth != gen.end(); ++truth)
480 {
481 if (truth->second != 0) // only signal event GenParticles
482 continue;
483
484 classification=-999;
485
486 auto particle = truth->first;
487
488 //determine eta region
489 Eta = std::abs(particle->momentum().pseudoRapidity());
490 if (Eta < m_cuts.maxEtaBarrel) Region = ETA_BARREL;
491 else if (Eta < m_cuts.maxEtaTransition) Region = ETA_TRANSITION;
492 else if (Eta < m_cuts.maxEtaEndcap) Region = ETA_ENDCAP;
493 else if ((Eta > m_cuts.minEtaFORWARD) && (Eta < m_cuts.maxEtaFORWARD)) Region = ETA_FORWARD;
494 else Region = ETA_OUTSIDE;
495 if ( !(Region==ETA_BARREL || Region==ETA_TRANSITION || Region==ETA_ENDCAP || Region==ETA_FORWARD) )
496 {
497 continue; // Only want to tally tracks that are within the detector volume.
498 }
499 tracks.m_counter[kTracks_gen_signal][TRACK_ALL][ETA_ALL] ++;
500 tracks.m_counter[kTracks_gen_signal][TRACK_ALL][Region]++;
501
502 //classify gen partilce as primary, secondary or truncated
503 classification = ClassifyParticle(particle,1.);
504
505 if(classification==TRACK_PRIMARY){
506 tracks.m_counter[kTracks_gen_signal][TRACK_PRIMARY][ETA_ALL] ++;
507 tracks.m_counter[kTracks_gen_signal][TRACK_PRIMARY][Region]++;
508 }
509 if(classification==TRACK_SECONDARY){
510 tracks.m_counter[kTracks_gen_signal][TRACK_SECONDARY][ETA_ALL] ++;
511 tracks.m_counter[kTracks_gen_signal][TRACK_SECONDARY][Region]++;
512 }
513
514 //see if gen track has at least 1 matching reco track with high enough probability
515 bool matched = false;
516 int nmatched = 0;
517
518 HepMcParticleLink hmpl2(particle,truth->second,HepMcParticleLink::IS_EVENTNUM); // FIXME truth->second is actually the position of the GenEvent in the McEventCollection!! See InDetRecStatisticsAlg::selectGenSignal(...) method (only client of TrackStatsHelper)
519 recoToTruthMap::iterator rttIter=rttMap.find(hmpl2);
520 if(rttIter != rttMap.end()){
521 for(imap = rttMap.lower_bound(hmpl2); imap !=rttMap.upper_bound(hmpl2); ++imap){
522 if(imap->second > m_cuts.matchTrackCut){
523 matched = true;
524 nmatched++;
525 }
526 }
527 }
528 if (matched) {
529 tracks.m_counter[kTracks_gen_signal][TRACK_MATCHED][Region]++;
530 tracks.m_counter[kTracks_gen_signal][TRACK_MATCHED][ETA_ALL]++;
531 if(nmatched > 1){
532 tracks.m_counter[kTracks_gen_signal][TRACK_MULTMATCH][Region] += nmatched-1;
533 tracks.m_counter[kTracks_gen_signal][TRACK_MULTMATCH][ETA_ALL]+= nmatched-1;
534 }
535 if(classification==TRACK_PRIMARY){
536 tracks.m_counter[kTracks_gen_signal][TRACK_MATCHED_PRIMARY][Region]++;
537 tracks.m_counter[kTracks_gen_signal][TRACK_MATCHED_PRIMARY][ETA_ALL]++;
538 if(nmatched > 1){
539 tracks.m_counter[kTracks_gen_signal][TRACK_MULTMATCH_PRIMARY][Region] += nmatched-1;
540 tracks.m_counter[kTracks_gen_signal][TRACK_MULTMATCH_PRIMARY][ETA_ALL] += nmatched-1;
541 }
542 }
543 if(classification==TRACK_SECONDARY){
544 tracks.m_counter[kTracks_gen_signal][TRACK_MATCHED_SECONDARY][Region]++;
545 tracks.m_counter[kTracks_gen_signal][TRACK_MATCHED_SECONDARY][ETA_ALL]++;
546 if(nmatched > 1){
547 tracks.m_counter[kTracks_gen_signal][TRACK_MULTMATCH_SECONDARY][Region] += nmatched-1;
548 tracks.m_counter[kTracks_gen_signal][TRACK_MULTMATCH_SECONDARY][ETA_ALL] += nmatched-1;
549 }
550 }
551 }
552 }
553
554 if (truth_missing) {
555 m_truthMissing = truth_missing;
556 }
557 m_tracks += tracks;
558 m_hits += hits;
559 m_trackSummarySum += trackSummarySum;
560 for (unsigned int i=0; i < Trk::TrackInfo::NumberOfTrackFitters; ++i) {
561 if (author_found[i]) {
562 m_author_found[i] = author_found[i];
563 }
564 }
565 {
566 std::lock_guard<std::mutex> lock(m_authorMutex);
567 m_recoInfo |= reco_info;
568 m_patternProperties |= pattern_properties;
569 }
570}
571
573 m_events = 0;
574 m_hits.reset();
575 m_tracks.reset();
576 m_trackSummarySum.reset();
577
578 for (int i=0; i<Trk::TrackInfo::NumberOfTrackFitters; i++) {
579 m_author_found[i] = false;
580 }
581 m_truthMissing=false;
582}
583
584
585void InDet::TrackStatHelper::print(MsgStream &out) const {
586
587 out << " Printing Statistics for " << key();
588
589 out << "TrackCollection \"" << key() << "\" " << "\n" ;
590
591 if (m_events > 0) {
592 out << "(TrackAuthors:";
593 std::vector<std::string> author_string;
594 author_string.reserve(Trk::TrackInfo::NumberOfTrackFitters);
595 for (int i=0; i<Trk::TrackInfo::NumberOfTrackFitters; ++i) {
596 author_string.push_back( Trk::TrackInfo(static_cast<Trk::TrackInfo::TrackFitter>(i),
598 m_patternProperties,
599 m_recoInfo).dumpInfo() );
600 }
601 for (int i=0; i<Trk::TrackInfo::NumberOfTrackFitters; i++){
602 if (m_author_found[i]){
603 out << " " << author_string[i];
604 }
605 }
606 out << " )" << std::endl
607 << "TrackTruthCollection \"" << Truthkey() << "\"" << std::endl;
608
609 if (m_truthMissing)
610 {
612 out << " WARNING: TrackTruth missing for part of this TrackCollection, no efficiencies or fakerates included.!" << std::endl;
613 else
614 out << " INFO: Intentionally no TrackTruth for this TrackCollection, no efficiencies or fakerates included." << std::endl;
615 }
616 if (!m_truthMissing)
617 {
618 out << " \t\t\t ................................tracks................................" << std::endl;
619 out << " \t\t\tn/event\tSEff.pr\tSEff.sd\tSLowP1\tSLowP2\tLowP1\tLowP2\tnoLink\tmultM\t" << std::endl;
620 out << " total "; printRegion1(out,ETA_ALL);
621 out << " in barrel "; printRegion1(out,ETA_BARREL);
622 out << " in trans. "; printRegion1(out,ETA_TRANSITION);
623 out << " in endcap "; printRegion1(out,ETA_ENDCAP);
624 out << " in forwa."; printRegion1(out,ETA_FORWARD);
625
626 }
627 else // no truth
628 {
629 out << "\t" << "tracks, n/event:" << std::endl;
630 out << "\t\t" << "total" << std::setiosflags(std::ios::fixed | std::ios::showpoint) <<
631 std::setw(7) << std::setprecision(2) << "\t" << m_tracks.m_counter[kTracks_rec][TRACK_ALL][ETA_ALL]/(float) m_events << std::endl << std::setprecision(-1);
632 out << "\t\t" << "in barrel" << std::setiosflags(std::ios::fixed | std::ios::showpoint) <<
633 std::setw(7) << std::setprecision(2) << "\t" << m_tracks.m_counter[kTracks_rec][TRACK_ALL][ETA_BARREL]/(float) m_events << std::endl << std::setprecision(-1);
634 out << "\t\t" << "in trans." << std::setiosflags(std::ios::fixed | std::ios::showpoint) <<
635 std::setw(7) << std::setprecision(2) << "\t" << m_tracks.m_counter[kTracks_rec][TRACK_ALL][ETA_TRANSITION]/(float) m_events << std::endl << std::setprecision(-1);
636 out << "\t\t" << "in endcap" << std::setiosflags(std::ios::fixed | std::ios::showpoint) <<
637 std::setw(7) << std::setprecision(2) << "\t" << m_tracks.m_counter[kTracks_rec][TRACK_ALL][ETA_ENDCAP]/(float) m_events << std::endl << std::setprecision(-1);
638 out << "\t\t" << "in forwa." << std::setiosflags(std::ios::fixed | std::ios::showpoint) <<
639 std::setw(7) << std::setprecision(2) << "\t" << m_tracks.m_counter[kTracks_rec][TRACK_ALL][ETA_FORWARD]/(float) m_events << std::endl << std::setprecision(-1);
640 }
641 out << " \t\t\t ....................................hits/track............................" << std::endl;
642 out << " \t\t\t\ttotal\tPIX1\tPIX2\tPIX3\tSCT1\tSCT2\tSCT3\tSCT4\tSCT5to9\tStraws" << std::endl;
643 out << " total "; printRegion2(out,ETA_ALL, (float) m_tracks.m_counter[kTracks_rec][TRACK_ALL][ETA_ALL]);
644 out << " in barrel "; printRegion2(out,ETA_BARREL, (float) m_tracks.m_counter[kTracks_rec][TRACK_ALL][ETA_BARREL]);
645 out << " in trans. "; printRegion2(out,ETA_TRANSITION, (float) m_tracks.m_counter[kTracks_rec][TRACK_ALL][ETA_TRANSITION] );
646 out << " in endcap "; printRegion2(out,ETA_ENDCAP, (float) m_tracks.m_counter[kTracks_rec][TRACK_ALL][ETA_ENDCAP] );
647 out << " in forwa. "; printRegion2(out,ETA_FORWARD, (float) m_tracks.m_counter[kTracks_rec][TRACK_ALL][ETA_FORWARD] );
648
649 }
650 else
651 out << ": NO EVENTS PROCESSED! " << std::endl;
652}
653
654void InDet::TrackStatHelper::printRegion1(MsgStream &out, enum eta_region region) const {
655 out << std::setiosflags(std::ios::fixed | std::ios::showpoint) << std::setw(7) << std::setprecision(2)
656 << "\t" << m_tracks.m_counter[kTracks_rec][TRACK_ALL][region]/(float) m_events;
657
658 if (m_tracks.m_counter[kTracks_gen_signal][TRACK_PRIMARY][region]) {
659 out << "\t" << std::setprecision(4)
660 << m_tracks.m_counter[kTracks_gen_signal][TRACK_MATCHED_PRIMARY][region]/ (float) m_tracks.m_counter[kTracks_gen_signal][TRACK_PRIMARY][region] // track efficiency wrt. gen primary in signal
661 << std::setprecision(2);
662 }
663 else {
664 out << "\t" << "n/a" ;
665 }
666 if (m_tracks.m_counter[kTracks_gen_signal][TRACK_SECONDARY][region]) {
667 out << "\t" << std::setprecision(4)
668 << m_tracks.m_counter[kTracks_gen_signal][TRACK_MATCHED_SECONDARY][region]/ (float) m_tracks.m_counter[kTracks_gen_signal][TRACK_SECONDARY][region] // track efficiency wrt. gen secondary in signal
669 << std::setprecision(2);
670 }
671 else {
672 out << "\t" << "n/a" ;
673 }
674
675 if (m_tracks.m_counter[kTracks_rec][TRACK_ALL_SIGNAL][region]) {
676 out << "\t" << 100*(m_tracks.m_counter[kTracks_rec][TRACK_LOWTRUTHPROB_SIGNAL][region]/ (float) m_tracks.m_counter[kTracks_rec][TRACK_ALL_SIGNAL][region]) << "%" // track fake rate for signal only
677 << "\t" << 100*(m_tracks.m_counter[kTracks_rec][TRACK_LOWTRUTHPROB2_SIGNAL][region]/ (float) m_tracks.m_counter[kTracks_rec][TRACK_ALL_SIGNAL][region]) << "%"; // track fake rate for signal only
678 }
679 else {
680 out << "\t" << "n/a\tn/a" ;
681 }
682
683 if (m_tracks.m_counter[kTracks_rec][TRACK_ALL][region]) {
684 out << "\t" << 100*(m_tracks.m_counter[kTracks_rec][TRACK_LOWTRUTHPROB][region]/ (float) m_tracks.m_counter[kTracks_rec][TRACK_ALL][region]) << "%" // track fake rate
685 << "\t" << 100*(m_tracks.m_counter[kTracks_rec][TRACK_LOWTRUTHPROB2][region]/ (float) m_tracks.m_counter[kTracks_rec][TRACK_ALL][region]) << "%" // track fake rate
686 << "\t" << 100*(m_tracks.m_counter[kTracks_rec][TRACK_NOHEPMCPARTICLELINK][region]/ (float) m_tracks.m_counter[kTracks_rec][TRACK_ALL][region]) << "%" // rate of tracks without HepMcParticleLink
687 << "\t" << 100*(m_tracks.m_counter[kTracks_gen][TRACK_MULTMATCH][region]/ (float) m_tracks.m_counter[kTracks_rec][TRACK_ALL][region]) << "%"; // rate of tracks matched to the same truth track
688
689 }
690 else{
691 out << "\tn/a\tn/a\tn/a";
692 }
693 out << std::endl << std::setprecision(-1);
694}
695
696void InDet::TrackStatHelper::printRegion2(MsgStream &out, enum eta_region region, float denominator) const {
697 out << std::setiosflags(std::ios::fixed | std::ios::showpoint) << std::setw(7) << std::setprecision(2);
698 if (denominator > 0) {
699 out << std::setiosflags(std::ios::fixed | std::ios::showpoint) << std::setw(4) << std::setprecision(1)
700 <<" \t\t" << float(m_hits.m_counter[kHits_rec][HIT_ALL ][region]/denominator) << std::setprecision(2)
701 << "\t" << float(m_hits.m_counter[kHits_rec][HIT_PIX1 ][region]/denominator)
702 << "\t" << float(m_hits.m_counter[kHits_rec][HIT_PIX2 ][region]/denominator)
703 << "\t" << float(m_hits.m_counter[kHits_rec][HIT_PIX3 ][region]/denominator)
704 << "\t" << float(m_hits.m_counter[kHits_rec][HIT_SCT1 ][region]/denominator)
705 << "\t" << float(m_hits.m_counter[kHits_rec][HIT_SCT2 ][region]/denominator)
706 << "\t" << float(m_hits.m_counter[kHits_rec][HIT_SCT3 ][region]/denominator)
707 << "\t" << float(m_hits.m_counter[kHits_rec][HIT_SCT4 ][region]/denominator)
708 << "\t" << float(m_hits.m_counter[kHits_rec][HIT_SCT5TO9][region]/denominator)
709 << "\t" << float(m_hits.m_counter[kHits_rec][HIT_TRT_ALL][region]/denominator)
710 << std::endl << std::setprecision(-1);
711 }
712 else {
713 out << " Unable to calculate: Denominator=0." << std::endl << std::setprecision(-1);
714 }
715}
716
718 enum track_types track_type ,
719 enum eta_region eta_region ) const
720{
721 long ex_denom = m_trackSummarySum.m_counter[kNTrackSummaryOK][track_type][eta_region][kNumberOfPixelHits]
722 + m_trackSummarySum.m_counter[kNTrackSummaryBAD][track_type][eta_region][kNumberOfPixelHits];
723 if (ex_denom > 0)
724 {
725 out << std::setiosflags(std::ios::fixed | std::ios::showpoint)
726 << std::setw(8) << track_types_string[track_type]
727 << std::setw(25) << key();
728
729 if (m_events)
730 out << std::setw(6) << std::setprecision(2) << m_tracks.m_counter[kTracks_rec][track_type][eta_region]/(float) m_events;
731 else
732 out << std::setw(6) << "n/a";
733 for (unsigned int stype_i=0; stype_i< kNSummaryTypes; ++stype_i) {
734 out << std::setw(std::max(6,static_cast<int>(strlen(s_summaryTypeName[stype_i])+1)))
735 << std::setprecision(2); printTrackSummaryAverage(out, track_type , eta_region, stype_i );
736 }
737 out << std::endl << std::setprecision(-1);
738 return true; // yes, printed output
739 }
740 else
741 return false; // no, didn't print output
742}
743
745 enum track_types track_type ,
746 enum eta_region eta_region,
747 int summary_type ) const
748{
749 assert( summary_type < kNSummaryTypes);
750 long denom = m_trackSummarySum.m_counter[kNTrackSummaryOK][track_type][eta_region][summary_type]
751 + m_trackSummarySum.m_counter[kNTrackSummaryBAD][track_type][eta_region][summary_type];
752
753 if (denom > 0) {
754 out << float ( m_trackSummarySum.m_counter[kTrackSummarySum][track_type][eta_region][summary_type]/ (float) denom);
755 }
756 else {
757 out << "n/a";
758 }
759}
760
761void InDet::TrackStatHelper::printSecondary(MsgStream &out) const {
762 out << "TrackCollection \"" << key() << "\" " << std::endl;
763 if (m_events > 0) {
764 out << "(TrackAuthors:";
765 std::vector<std::string> author_string;
766 author_string.reserve(Trk::TrackInfo::NumberOfTrackFitters);
767 for (int i=0; i<Trk::TrackInfo::NumberOfTrackFitters; i++) {
768 author_string.push_back( Trk::TrackInfo(static_cast<Trk::TrackInfo::TrackFitter>(i),
770 m_patternProperties,
771 m_recoInfo).dumpInfo() );
772 }
773 for (int i=0; i<Trk::TrackInfo::NumberOfTrackFitters; i++){
774 if (m_author_found[i]){
775 out << " " << author_string[i];
776 }
777 }
778 out << " )" << std::endl
779 << "TrackTruthCollection \"" << Truthkey() << "\"" << std::endl;
780
781 if (m_truthMissing)
782 {
784 out << " WARNING: TrackTruth missing for part of this TrackCollection --> No secondaries information printed!" << std::endl;
785 else
786 out << " INFO: Intentionally no TrackTruth for this TrackCollection. (No secondaries information printed.)" << std::endl;
787 }
788 if (!m_truthMissing)
789 {
790 out << " \t\t\t\t ......................truth mached tracks statistics....................." << std::endl;
791 out << " \t\t\t\tn/event\teff.\ttotal\tPIX1\tPIX2\tPIX3\tSCT1\tSCT2\tSCT3\tSCT4\tSCT5to9\tStraws" << std::endl;
792
793
794 out << " total secondaries ";
795 printRegionSecondary(out, ETA_ALL, (float) (m_tracks.m_counter[kTracks_gen][TRACK_MATCHED_SECONDARY][ETA_ALL] +m_tracks.m_counter[kTracks_gen][TRACK_MULTMATCH_SECONDARY][ETA_ALL]) );
796 out << " secondaries in barrel ";
797 printRegionSecondary(out, ETA_BARREL, (float) (m_tracks.m_counter[kTracks_gen][TRACK_MATCHED_SECONDARY][ETA_BARREL] +m_tracks.m_counter[kTracks_gen][TRACK_MULTMATCH_SECONDARY][ETA_BARREL]) );
798 out << " secondaries in trans. ";
799 printRegionSecondary(out, ETA_TRANSITION, (float) (m_tracks.m_counter[kTracks_gen][TRACK_MATCHED_SECONDARY][ETA_TRANSITION]+m_tracks.m_counter[kTracks_gen][TRACK_MULTMATCH_SECONDARY][ETA_TRANSITION]) );
800 out << " secondaries in endcap ";
801 printRegionSecondary(out, ETA_ENDCAP, (float) (m_tracks.m_counter[kTracks_gen][TRACK_MATCHED_SECONDARY][ETA_ENDCAP] +m_tracks.m_counter[kTracks_gen][TRACK_MULTMATCH_SECONDARY][ETA_ENDCAP]) );
802 out << " secondaries in forwa. ";
803 printRegionSecondary(out, ETA_FORWARD, (float) (m_tracks.m_counter[kTracks_gen][TRACK_MATCHED_SECONDARY][ETA_FORWARD] +m_tracks.m_counter[kTracks_gen][TRACK_MULTMATCH_SECONDARY][ETA_FORWARD]) );
804 }
805 }
806 else out << ": NO EVENTS PROCESSED! " << std::endl;
807}
808
810 enum eta_region region,
811 float denominator) const {
812 out << std::setiosflags(std::ios::fixed | std::ios::showpoint) << std::setw(7) << std::setprecision(2)
813 << "\t" << m_tracks.m_counter[kTracks_gen][TRACK_MATCHED_SECONDARY][region]/(float) m_events; // tracks / event
814 if (m_tracks.m_counter[kTracks_gen][TRACK_SECONDARY][region]) {
815 out << "\t" << m_tracks.m_counter[kTracks_gen][TRACK_MATCHED_SECONDARY][region]/ (float) m_tracks.m_counter[kTracks_gen][TRACK_SECONDARY][region]; // track efficiency
816 }
817 else {
818 out << "\t" << "n/a" ;
819 }
820
821 if (denominator > 0) {
822 out << std::setiosflags(std::ios::fixed | std::ios::showpoint) << std::setw(4) << std::setprecision(1)
823 <<" " << float(m_hits.m_counter[kHits_sec][HIT_ALL ][region]/denominator) << std::setprecision(2)
824 << "\t" << float(m_hits.m_counter[kHits_sec][HIT_PIX1 ][region]/denominator)
825 << "\t" << float(m_hits.m_counter[kHits_sec][HIT_PIX2 ][region]/denominator)
826 << "\t" << float(m_hits.m_counter[kHits_sec][HIT_PIX3 ][region]/denominator)
827 << "\t" << float(m_hits.m_counter[kHits_sec][HIT_SCT1 ][region]/denominator)
828 << "\t" << float(m_hits.m_counter[kHits_sec][HIT_SCT2 ][region]/denominator)
829 << "\t" << float(m_hits.m_counter[kHits_sec][HIT_SCT3 ][region]/denominator)
830 << "\t" << float(m_hits.m_counter[kHits_sec][HIT_SCT4 ][region]/denominator)
831 << "\t" << float(m_hits.m_counter[kHits_sec][HIT_SCT5TO9][region]/denominator)
832 << "\t" << float(m_hits.m_counter[kHits_sec][HIT_TRT_ALL][region]/denominator)
833 << std::endl << std::setprecision(-1);
834 }
835 else {
836 out << " Unable to calculate: Denominator=0." << std::endl << std::setprecision(-1);
837 }
838}
839
840
841
843 bool passed = false;
844 if(para->pT() > m_cuts.minPt && std::abs(para->eta()) < m_cuts.maxEtaFORWARD)passed = true;
845
846
847 return passed;
848
849}
850
851int InDet::TrackStatHelper::ClassifyParticle( const HepMC::ConstGenParticlePtr& particle, const double prob) const {
852
853 int partClass=-999;
854
855 // ME: only classify matched tracks
856 if (prob <= m_cuts.matchTrackCut) return partClass;
857
858 //classify as primary, secondary or truncated
859 bool primary=false;
860 bool secondary=false;
861 bool truncated=false;
862
863 if (particle->production_vertex()){
864
865 // primary vertex inside innermost layer?
866 HepGeom::Point3D<double> startVertex(particle->production_vertex()->position().x(),
867 particle->production_vertex()->position().y(),
868 particle->production_vertex()->position().z());
869 if ( std::abs(startVertex.perp()) < m_cuts.maxRStartPrimary
870 && std::abs(startVertex.z()) < m_cuts.maxZStartPrimary) {
871 if (particle->end_vertex() == nullptr) {
872 primary=true;
873 }
874 else {
875 HepGeom::Point3D<double> endVertex(particle->end_vertex()->position().x(),
876 particle->end_vertex()->position().y(),
877 particle->end_vertex()->position().z());
878 if ( endVertex.perp() > m_cuts.minREndPrimary
879 || std::abs(startVertex.z()) > m_cuts.minZEndPrimary){
880 primary=true;
881 }
882 else {
883 truncated = true;
884 }
885 }
886 }
887
888 else if ( startVertex.perp() < m_cuts.maxRStartSecondary && std::abs(startVertex.z()) < m_cuts.maxZStartSecondary) {
889 if (particle->end_vertex() == nullptr) {
890 secondary=true;
891 }
892 else {
893 HepGeom::Point3D<double> endVertex(particle->end_vertex()->position().x(),
894 particle->end_vertex()->position().y(),
895 particle->end_vertex()->position().z());
896
897 if (endVertex.perp() > m_cuts.minREndSecondary
898 || std::abs(endVertex.z()) > m_cuts.minZEndSecondary) {
899 secondary=true;
900 }
901 }
902 }
903 }//end classification
904
905 if(secondary)partClass=TRACK_SECONDARY;
906 if(primary)partClass=TRACK_PRIMARY;
907 if(truncated)partClass=TRACK_TRUNCATED;
908
909 return partClass;
910}
This class provides an interface to generate or decode an identifier for the upper levels of the dete...
virtual void lock()=0
Interface to allow an object to lock itself when made const in SG.
bool passed(DecisionID id, const DecisionIDContainer &)
checks if required decision ID is in the set of IDs in the container
This is an Identifier helper class for the Pixel subdetector.
This is an Identifier helper class for the SCT subdetector.
DataVector< Trk::Track > TrackCollection
This typedef represents a collection of Trk::Track objects.
static const std::string track_types_string[InDet::N_TRACKTYPES]
This class provides an interface to generate or decode an identifier for the upper levels of the dete...
bool is_sct(Identifier id) const
bool is_pixel(Identifier id) const
bool is_trt(Identifier id) const
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.
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
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.
@ 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 ...
This is an Identifier helper class for the Pixel subdetector.
Definition PixelID.h:69
int layer_disk(const Identifier &id) const
Definition PixelID.h:602
This is an Identifier helper class for the SCT subdetector.
Definition SCT_ID.h:68
int layer_disk(const Identifier &id) const
Definition SCT_ID.h:687
MC particle associated with a reco track + the quality of match.
Definition TrackTruth.h:14
float probability() const
Definition TrackTruth.h:28
const HepMcParticleLink & particleLink() const
Definition TrackTruth.h:26
Base class for all CompetingRIOsOnTack implementations, extends the common MeasurementBase.
Interface for condensing Trk::Track properties and associated hits to a (non-fittable) foot print,...
virtual std::unique_ptr< Trk::TrackSummary > summary(const EventContext &ctx, const Track &track) const =0
Start from a copy of the existing input track summary if there, otherwise start from a new one.
This class is the pure abstract base class for all fittable tracking measurements.
virtual bool type(MeasurementBaseType::Type type) const =0
Interface method checking the type.
double eta() const
Access method for pseudorapidity - from momentum.
double pT() const
Access method for transverse momentum.
Class to handle RIO On Tracks ROT) for InDet and Muons, it inherits from the common MeasurementBase.
Definition RIO_OnTrack.h:70
Identifier identify() const
return the identifier -extends MeasurementBase
Small utility to cast MeasurementBase (either in a vector or single) into RIO_OnTrack.
static void extract(std::vector< const RIO_OnTrack * > &rots, const std::vector< const MeasurementBase * > &measurements)
Contains information about the 'fitter' of this track.
TrackFitter
enums to identify who created this track and what propertis does it have.
represents the track state (measurement, material, fit parameters and quality) at a surface.
@ Measurement
This is a measurement, and will at least contain a Trk::MeasurementBase.
A summary of the information contained by a track.
const GenParticle * ConstGenParticlePtr
Definition GenParticle.h:38
Ensure that the ATLAS eigen extensions are properly loaded.
ParametersBase< TrackParametersDim, Charged > TrackParameters
SummaryType
enumerates the different types of information stored in Summary.
@ numberOfGangedPixels
number of Ganged Pixels flagged as fakes
@ numberOfPixelHits
number of pixel layers on track with absence of hits
@ numberOfTRTHighThresholdOutliers
number of dead TRT straws crossed
@ numberOfInnermostPixelLayerHits
these are the hits in the 1st pixel layer
@ numberOfTRTHighThresholdHits
total number of TRT hits which pass the high threshold
@ numberOfSCTHoles
number of Holes in both sides of a SCT module
@ numberOfPixelHoles
number of pixels which have a ganged ambiguity.
@ numberOfOutliersOnTrack
100 times the standard deviation of the chi2 from the surfaces
STL namespace.
T_Int m_counter[N_Categories][N_Types][N_Regions]