CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Private Attributes
TrackingParticleNumberOfLayers Class Reference

#include <TrackingParticleNumberOfLayers.h>

Public Types

enum  { nTrackerLayers = 0, nPixelLayers = 1, nStripMonoAndStereoLayers = 2 }
 

Public Member Functions

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
 
 TrackingParticleNumberOfLayers (const edm::Event &iEvent, const std::vector< edm::EDGetTokenT< std::vector< PSimHit >>> &simHitTokens)
 

Private Attributes

std::vector< std::pair< unsigned int, const PSimHit * > > trackIdToHitPtr_
 

Detailed Description

This class calculates the number of tracker layers, pixel layers, and strip mono+stereo layers "crossed" by TrackingParticle.

The numbers of pixel and strip mono+stereo layers are not available from TrackingParticle itself, so they are calculated here in a standalone way in order to not to modify the TP dataformat (for now). The number of tracker layers is available in TP, but its calculation in TrackingTruthAccumulator gives wrong results (too many layers) for loopers, so also it is calculated here on the same go.

The PSimHits are needed for the calculation, so, in practice, in events with pileup the numbers of layers can be calculated only for TPs from the signal event (i.e. not for pileup TPs). Fortunately this is exactly what is sufficient for MultiTrackValidator.

Eventually we should move to use HitPattern as in reco::TrackBase (more information in a compact format), and consider adding it to the TrackingParticle itself.

In principle we could utilize the TP->SimHit map produced in SimHitTPAssociationProducer instead of doing the association here again, but

Definition at line 46 of file TrackingParticleNumberOfLayers.h.

Member Enumeration Documentation

◆ anonymous enum

anonymous enum

Constructor & Destructor Documentation

◆ TrackingParticleNumberOfLayers()

TrackingParticleNumberOfLayers::TrackingParticleNumberOfLayers ( const edm::Event iEvent,
const std::vector< edm::EDGetTokenT< std::vector< PSimHit >>> &  simHitTokens 
)

Definition at line 28 of file TrackingParticleNumberOfLayers.cc.

References iEvent, rpcPointValidation_cfi::simHit, SiPixelPhase1TrackingParticleV_cfi::simHitToken, and trackIdToHitPtr_.

29  {
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 }
int iEvent
Definition: GenABIO.cc:224
std::vector< std::pair< unsigned int, const PSimHit * > > trackIdToHitPtr_

Member Function Documentation

◆ calculate()

std::tuple< std::unique_ptr< edm::ValueMap< unsigned int > >, std::unique_ptr< edm::ValueMap< unsigned int > >, std::unique_ptr< edm::ValueMap< unsigned int > > > TrackingParticleNumberOfLayers::calculate ( const edm::Handle< TrackingParticleCollection > &  tps,
const TrackerTopology tTopo 
) const

Definition at line 46 of file TrackingParticleNumberOfLayers.cc.

References ALPAKA_ACCELERATOR_NAMESPACE::brokenline::constexpr(), trigObjTnPSource_cfi::filler, mps_fire::i, fastTrackerRecHitType::isPixel(), dqmiolumiharvest::j, nano_mu_digi_cff::layer, TrackerTopology::layer(), eostools::move(), MuonTCMETValueMapProducer_cff::nLayers, nPixelLayers, nStripMonoAndStereoLayers, PbPb_ZMuSkimMuonDPG_cff::particleType, EgammaValidation_cff::pdgId, PixelSubdetector::PixelBarrel, PixelSubdetector::PixelEndcap, FastTimerService_cff::range, rpcPointValidation_cfi::simHit, cscDigiValidation_cfi::simTrack, DetId::subdetId(), StripSubdetector::TEC, TrackerTopology::tecIsRPhi(), TrackerTopology::tecIsStereo(), StripSubdetector::TIB, TrackerTopology::tibIsRPhi(), TrackerTopology::tibIsStereo(), StripSubdetector::TID, TrackerTopology::tidIsRPhi(), TrackerTopology::tidIsStereo(), StripSubdetector::TOB, TrackerTopology::tobIsRPhi(), TrackerTopology::tobIsStereo(), cmsswSequenceInfo::tp, and trackIdToHitPtr_.

47  {
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
bool tibIsStereo(const DetId &id) const
unsigned int layer(const DetId &id) const
bool tobIsStereo(const DetId &id) const
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: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
bool tecIsStereo(const DetId &id) const
bool tidIsStereo(const DetId &id) const
bool isPixel(HitType hitType)
Monte Carlo truth information used for tracking validation.
std::vector< TrackingParticle > TrackingParticleCollection
static constexpr auto TID
def move(src, dest)
Definition: eostools.py:511

Member Data Documentation

◆ trackIdToHitPtr_

std::vector<std::pair<unsigned int, const PSimHit *> > TrackingParticleNumberOfLayers::trackIdToHitPtr_
private

Definition at line 59 of file TrackingParticleNumberOfLayers.h.

Referenced by calculate(), and TrackingParticleNumberOfLayers().