CMS 3D CMS Logo

TrackingParticleNumberOfLayers.cc
Go to the documentation of this file.
2 
6 
7 namespace {
8  bool trackIdHitPairLess(const std::pair<unsigned int, const PSimHit *> &a,
9  const std::pair<unsigned int, const PSimHit *> &b) {
10  return a.first < b.first;
11  }
12 
13  bool trackIdHitPairLessSort(const std::pair<unsigned int, const PSimHit *> &a,
14  const std::pair<unsigned int, const PSimHit *> &b) {
15  if (a.first == b.first) {
16  const auto atof = edm::isFinite(a.second->timeOfFlight())
17  ? a.second->timeOfFlight()
18  : std::numeric_limits<decltype(a.second->timeOfFlight())>::max();
19  const auto btof = edm::isFinite(b.second->timeOfFlight())
20  ? b.second->timeOfFlight()
21  : std::numeric_limits<decltype(b.second->timeOfFlight())>::max();
22  return atof < btof;
23  }
24  return a.first < b.first;
25  }
26 } // namespace
27 
29  const edm::Event &iEvent, const std::vector<edm::EDGetTokenT<std::vector<PSimHit>>> &simHitTokens) {
30  // A multimap linking SimTrack::trackId() to a pointer to PSimHit
31  // Similar to TrackingTruthAccumulator
32  for (const auto &simHitToken : simHitTokens) {
34  iEvent.getByToken(simHitToken, hsimhits);
35  trackIdToHitPtr_.reserve(trackIdToHitPtr_.size() + hsimhits->size());
36  for (const auto &simHit : *hsimhits) {
37  trackIdToHitPtr_.emplace_back(simHit.trackId(), &simHit);
38  }
39  }
40  std::stable_sort(trackIdToHitPtr_.begin(), trackIdToHitPtr_.end(), trackIdHitPairLessSort);
41 }
42 
43 std::tuple<std::unique_ptr<edm::ValueMap<unsigned int>>,
44  std::unique_ptr<edm::ValueMap<unsigned int>>,
45  std::unique_ptr<edm::ValueMap<unsigned int>>>
47  const TrackerTopology &tTopo) const {
48  const TrackingParticleCollection &tps = *htps;
49  std::vector<unsigned int> valuesLayers(tps.size(), 0);
50  std::vector<unsigned int> valuesPixelLayers(tps.size(), 0);
51  std::vector<unsigned int> valuesStripMonoAndStereoLayers(tps.size(), 0);
52  for (size_t iTP = 0; iTP < tps.size(); ++iTP) {
53  const TrackingParticle &tp = tps[iTP];
54  const auto pdgId = tp.pdgId();
55 
56  // I would prefer a better way...
57  constexpr unsigned int maxSubdet = static_cast<unsigned>(StripSubdetector::TEC) + 1;
58  constexpr unsigned int maxLayer = 0xF + 1; // as in HitPattern.h
59  bool hasHit[maxSubdet][maxLayer];
60  bool hasPixel[maxSubdet][maxLayer];
61  bool hasMono[maxSubdet][maxLayer];
62  bool hasStereo[maxSubdet][maxLayer];
63  memset(hasHit, 0, sizeof(hasHit));
64  memset(hasPixel, 0, sizeof(hasPixel));
65  memset(hasMono, 0, sizeof(hasMono));
66  memset(hasStereo, 0, sizeof(hasStereo));
67 
68  for (const SimTrack &simTrack : tp.g4Tracks()) {
69  // Logic is from TrackingTruthAccumulator
70  auto range = std::equal_range(trackIdToHitPtr_.begin(),
71  trackIdToHitPtr_.end(),
72  std::pair<unsigned int, const PSimHit *>(simTrack.trackId(), nullptr),
73  trackIdHitPairLess);
74  if (range.first == range.second)
75  continue;
76 
77  auto iHitPtr = range.first;
78  for (; iHitPtr != range.second; ++iHitPtr) {
79  // prevent simhits with particleType != pdgId from being the "first hit"
80  if (iHitPtr->second->particleType() == pdgId)
81  break;
82  }
83  if (iHitPtr == range.second) // no simhits with particleType == pdgId
84  continue;
85  int processType = iHitPtr->second->processType();
86  int particleType = iHitPtr->second->particleType();
87 
88  for (; iHitPtr != range.second; ++iHitPtr) {
89  const PSimHit &simHit = *(iHitPtr->second);
90  if (simHit.eventId() != tp.eventId())
91  continue;
92  DetId newDetector = DetId(simHit.detUnitId());
93 
94  // Check for delta and interaction products discards
95  if (processType == simHit.processType() && particleType == simHit.particleType() &&
96  pdgId == simHit.particleType()) {
97  // The logic of this piece follows HitPattern
98  bool isPixel = false;
99  bool isStripStereo = false;
100  bool isStripMono = false;
101 
102  switch (newDetector.subdetId()) {
105  isPixel = true;
106  break;
108  isStripMono = tTopo.tibIsRPhi(newDetector);
109  isStripStereo = tTopo.tibIsStereo(newDetector);
110  break;
112  isStripMono = tTopo.tidIsRPhi(newDetector);
113  isStripStereo = tTopo.tidIsStereo(newDetector);
114  break;
116  isStripMono = tTopo.tobIsRPhi(newDetector);
117  isStripStereo = tTopo.tobIsStereo(newDetector);
118  break;
120  isStripMono = tTopo.tecIsRPhi(newDetector);
121  isStripStereo = tTopo.tecIsStereo(newDetector);
122  break;
123  }
124 
125  const auto subdet = newDetector.subdetId();
126  const auto layer = tTopo.layer(newDetector);
127 
128  hasHit[subdet][layer] = true;
129  if (isPixel)
130  hasPixel[subdet][layer] = isPixel;
131  else if (isStripMono)
132  hasMono[subdet][layer] = isStripMono;
133  else if (isStripStereo)
134  hasStereo[subdet][layer] = isStripStereo;
135  }
136  }
137  }
138 
139  unsigned int nLayers = 0;
140  unsigned int nPixelLayers = 0;
141  unsigned int nStripMonoAndStereoLayers = 0;
142  for (unsigned int i = 0; i < maxSubdet; ++i) {
143  for (unsigned int j = 0; j < maxLayer; ++j) {
144  nLayers += hasHit[i][j];
145  nPixelLayers += hasPixel[i][j];
146  nStripMonoAndStereoLayers += (hasMono[i][j] && hasStereo[i][j]);
147  }
148  }
149 
150  valuesLayers[iTP] = nLayers;
151  valuesPixelLayers[iTP] = nPixelLayers;
152  valuesStripMonoAndStereoLayers[iTP] = nStripMonoAndStereoLayers;
153  }
154 
155  auto ret0 = std::make_unique<edm::ValueMap<unsigned int>>();
156  {
158  filler.insert(htps, valuesLayers.begin(), valuesLayers.end());
159  filler.fill();
160  }
161  auto ret1 = std::make_unique<edm::ValueMap<unsigned int>>();
162  {
164  filler.insert(htps, valuesPixelLayers.begin(), valuesPixelLayers.end());
165  filler.fill();
166  }
167  auto ret2 = std::make_unique<edm::ValueMap<unsigned int>>();
168  {
170  filler.insert(htps, valuesStripMonoAndStereoLayers.begin(), valuesStripMonoAndStereoLayers.end());
171  filler.fill();
172  }
173 
174  return std::make_tuple(std::move(ret0), std::move(ret1), std::move(ret2));
175 }
static constexpr auto TEC
TrackingParticleNumberOfLayers(const edm::Event &iEvent, const std::vector< edm::EDGetTokenT< std::vector< PSimHit >>> &simHitTokens)
constexpr bool isFinite(T x)
bool tibIsStereo(const DetId &id) const
unsigned int layer(const DetId &id) const
bool tobIsStereo(const DetId &id) const
constexpr std::array< uint8_t, layerIndexSize< TrackerTraits > > layer
bool tibIsRPhi(const DetId &id) const
int iEvent
Definition: GenABIO.cc:224
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
static constexpr auto TOB
std::vector< std::pair< unsigned int, const PSimHit * > > trackIdToHitPtr_
bool tobIsRPhi(const DetId &id) const
bool tidIsRPhi(const DetId &id) const
bool tecIsRPhi(const DetId &id) const
Definition: DetId.h:17
static constexpr auto TIB
double b
Definition: hdecay.h:118
bool tecIsStereo(const DetId &id) const
bool tidIsStereo(const DetId &id) const
double a
Definition: hdecay.h:119
bool isPixel(HitType hitType)
Monte Carlo truth information used for tracking validation.
std::vector< TrackingParticle > TrackingParticleCollection
static constexpr auto TID
std::tuple< std::unique_ptr< edm::ValueMap< unsigned int > >, std::unique_ptr< edm::ValueMap< unsigned int > >, std::unique_ptr< edm::ValueMap< unsigned int > > > calculate(const edm::Handle< TrackingParticleCollection > &tps, const TrackerTopology &tTopo) const
def move(src, dest)
Definition: eostools.py:511