CMS 3D CMS Logo

SiPixelPhase1TrackingParticleV.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: SiPixelPhase1TrackingParticleV
4 // Class: SiPixelPhase1TrackingParticleV
5 //
6 
7 // Original Author: Marcel Schneider
8 // Additional Authors: Alexander Morton - modifying code for validation use
9 
12 
15 
18 
21 
22 class TrackAssociatorByHits;
23 
24 namespace {
25  bool trackIdHitPairLess(const std::pair<unsigned int, const PSimHit*>& a, const std::pair<unsigned int, const PSimHit*>& b) {
26  return a.first < b.first;
27  }
28 
29  bool trackIdHitPairLessSort(const std::pair<unsigned int, const PSimHit*>& a, const std::pair<unsigned int, const PSimHit*>& b) {
30  if(a.first == b.first) {
31  const auto atof = edm::isFinite(a.second->timeOfFlight()) ? a.second->timeOfFlight() : std::numeric_limits<decltype(a.second->timeOfFlight())>::max();
32  const auto btof = edm::isFinite(b.second->timeOfFlight()) ? b.second->timeOfFlight() : std::numeric_limits<decltype(b.second->timeOfFlight())>::max();
33  return atof < btof;
34  }
35  return a.first < b.first;
36  }
37 }
38 
39 
41  SiPixelPhase1Base(iConfig),
42  vec_TrackingParticle_Token_( consumes<TrackingParticleCollection>( iConfig.getParameter<edm::InputTag>( "src" ) ) )
43 {
44  for(const auto& tag: iConfig.getParameter<std::vector<edm::InputTag>>("simHitToken")) {
45  simHitTokens_.push_back(consumes<std::vector<PSimHit>>(tag));
46  }
47 }
48 
50 
51  edm::Handle<TrackingParticleCollection> TruthTrackContainer;
52  iEvent.getByToken( vec_TrackingParticle_Token_, TruthTrackContainer );
53  const TrackingParticleCollection *tPC = TruthTrackContainer.product();
54 
55  // A multimap linking SimTrack::trackId() to a pointer to PSimHit
56  // Similar to TrackingTruthAccumulator
57  for(const auto& simHitToken: simHitTokens_) {
59  iEvent.getByToken(simHitToken, hsimhits);
60  trackIdToHitPtr_.reserve(trackIdToHitPtr_.size()+hsimhits->size());
61  for(const auto& simHit: *hsimhits) {
62  trackIdToHitPtr_.emplace_back(simHit.trackId(), &simHit);
63  }
64  }
65  std::stable_sort(trackIdToHitPtr_.begin(), trackIdToHitPtr_.end(), trackIdHitPairLessSort);
66 
67 
68  // Loop over TrackingParticle's
69  for (TrackingParticleCollection::const_iterator t = tPC -> begin(); t != tPC -> end(); ++t) {
70 
71  // histo manager requires a det ID, use first tracker hit
72 
73  bool isBpixtrack = false, isFpixtrack = false;
74  DetId id;
75 
76  for(const SimTrack& simTrack: t->g4Tracks()) {
77  // Logic is from TrackingTruthAccumulator
78  auto range = std::equal_range(trackIdToHitPtr_.begin(), trackIdToHitPtr_.end(), std::pair<unsigned int, const PSimHit *>(simTrack.trackId(), nullptr), trackIdHitPairLess);
79  if(range.first == range.second) continue;
80 
81  auto iHitPtr = range.first;
82  for(; iHitPtr != range.second; ++iHitPtr) {
83  const PSimHit& simHit = *(iHitPtr->second);
84  if(simHit.eventId() != t->eventId())
85  continue;
86  id = DetId( simHit.detUnitId() );
87 
88  // check we are in pixel
89  uint32_t subdetid = (id.subdetId());
90  if (subdetid == PixelSubdetector::PixelBarrel) isBpixtrack = true;
91  if (subdetid == PixelSubdetector::PixelEndcap) isFpixtrack = true;
92  if (subdetid != PixelSubdetector::PixelBarrel && subdetid != PixelSubdetector::PixelEndcap) continue;
93  }
94  }
95 
96  if ( isBpixtrack || isFpixtrack ) {
97  histo[MASS].fill(t->mass(), id, &iEvent);
98  histo[CHARGE].fill(t->charge(), id, &iEvent);
99  histo[ID].fill(t->pdgId(), id, &iEvent);
100  histo[NHITS].fill(t->numberOfTrackerHits(), id, &iEvent);
101  histo[MATCHED].fill(t->numberOfTrackerLayers(), id, &iEvent);
102  histo[PT].fill(sqrt(t->momentum().perp2()), id, &iEvent);
103  histo[PHI].fill(t->momentum().Phi(), id, &iEvent);
104  histo[ETA].fill(t->momentum().eta(), id, &iEvent);
105  histo[VTX].fill(t->vx(), id, &iEvent);
106  histo[VTY].fill(t->vy(), id, &iEvent);
107  histo[VYZ].fill(t->vz(), id, &iEvent);
108  histo[TIP].fill(sqrt(t->vertex().perp2()), id, &iEvent);
109  histo[LIP].fill(t->vz(), id, &iEvent);
110  }
111  }
112 }
113 
115 
T getParameter(std::string const &) const
std::vector< TrackingParticle > TrackingParticleCollection
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
bool isFinite(T x)
int iEvent
Definition: GenABIO.cc:230
simTrack
per collection params
T sqrt(T t)
Definition: SSEVec.h:18
EncodedEventId eventId() const
Definition: PSimHit.h:105
#define end
Definition: vmac.h:37
std::vector< std::pair< unsigned int, const PSimHit * > > trackIdToHitPtr_
edm::EDGetTokenT< TrackingParticleCollection > vec_TrackingParticle_Token_
Definition: DetId.h:18
T const * product() const
Definition: Handle.h:81
double b
Definition: hdecay.h:120
void analyze(const edm::Event &, const edm::EventSetup &)
std::vector< HistogramManager > histo
#define begin
Definition: vmac.h:30
HLT enums.
double a
Definition: hdecay.h:121
SiPixelPhase1TrackingParticleV(const edm::ParameterSet &conf)
std::vector< edm::EDGetTokenT< std::vector< PSimHit > > > simHitTokens_
unsigned int detUnitId() const
Definition: PSimHit.h:93