CMS 3D CMS Logo

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