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 edm::EventSetup &iSetup) 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 44 of file TrackingParticleNumberOfLayers.h.

Member Enumeration Documentation

anonymous enum

Constructor & Destructor Documentation

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

Definition at line 24 of file TrackingParticleNumberOfLayers.cc.

References edm::Event::getByToken(), rpcPointValidation_cfi::simHit, SiPixelPhase1TrackingParticleV_cfi::simHitToken, and trackIdToHitPtr_.

24  {
25 
26  // A multimap linking SimTrack::trackId() to a pointer to PSimHit
27  // Similar to TrackingTruthAccumulator
28  for(const auto& simHitToken: simHitTokens) {
30  iEvent.getByToken(simHitToken, hsimhits);
31  trackIdToHitPtr_.reserve(trackIdToHitPtr_.size()+hsimhits->size());
32  for(const auto& simHit: *hsimhits) {
33  trackIdToHitPtr_.emplace_back(simHit.trackId(), &simHit);
34  }
35  }
36  std::stable_sort(trackIdToHitPtr_.begin(), trackIdToHitPtr_.end(), trackIdHitPairLessSort);
37 }
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
std::vector< std::pair< unsigned int, const PSimHit * > > trackIdToHitPtr_

Member Function Documentation

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 edm::EventSetup iSetup 
) const

Definition at line 42 of file TrackingParticleNumberOfLayers.cc.

References constexpr, PSimHit::detUnitId(), TrackingParticle::eventId(), PSimHit::eventId(), edm::helper::Filler< Map >::fill(), objects.autophobj::filler, TrackingParticle::g4Tracks(), edm::EventSetup::get(), mps_fire::i, edm::helper::Filler< Map >::insert(), fastTrackerRecHitType::isPixel(), TrackerTopology::layer(), eostools::move(), MuonTCMETValueMapProducer_cff::nLayers, nPixelLayers, nStripMonoAndStereoLayers, objects.autophobj::particleType, PSimHit::particleType(), common_cff::pdgId, TrackingParticle::pdgId(), PixelSubdetector::PixelBarrel, PixelSubdetector::PixelEndcap, PSimHit::processType(), rpcPointValidation_cfi::simHit, simTrackMatching_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(), and trackIdToHitPtr_.

Referenced by TrackingParticleNumberOfLayersProducer::produce().

42  {
43  //Retrieve tracker topology from geometry
45  iSetup.get<TrackerTopologyRcd>().get(tTopoHandle);
46  const TrackerTopology& tTopo = *tTopoHandle;
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(), trackIdToHitPtr_.end(), std::pair<unsigned int, const PSimHit *>(simTrack.trackId(), nullptr), trackIdHitPairLess);
71  if(range.first == range.second) continue;
72 
73  auto iHitPtr = range.first;
74  for(; iHitPtr != range.second; ++iHitPtr) {
75  // prevent simhits with particleType != pdgId from being the "first hit"
76  if(iHitPtr->second->particleType() == pdgId)
77  break;
78  }
79  if(iHitPtr == range.second) // no simhits with particleType == pdgId
80  continue;
81  int processType = iHitPtr->second->processType();
82  int particleType = iHitPtr->second->particleType();
83 
84  for(; iHitPtr != range.second; ++iHitPtr) {
85  const PSimHit& simHit = *(iHitPtr->second);
86  if(simHit.eventId() != tp.eventId())
87  continue;
88  DetId newDetector = DetId( simHit.detUnitId() );
89 
90  // Check for delta and interaction products discards
91  if( processType == simHit.processType() && particleType == simHit.particleType() && pdgId == simHit.particleType() ) {
92  // The logic of this piece follows HitPattern
93  bool isPixel = false;
94  bool isStripStereo = false;
95  bool isStripMono = false;
96 
97  switch(newDetector.subdetId()) {
100  isPixel = true;
101  break;
103  isStripMono = tTopo.tibIsRPhi(newDetector);
104  isStripStereo = tTopo.tibIsStereo(newDetector);
105  break;
107  isStripMono = tTopo.tidIsRPhi(newDetector);
108  isStripStereo = tTopo.tidIsStereo(newDetector);
109  break;
111  isStripMono = tTopo.tobIsRPhi(newDetector);
112  isStripStereo = tTopo.tobIsStereo(newDetector);
113  break;
115  isStripMono = tTopo.tecIsRPhi(newDetector);
116  isStripStereo = tTopo.tecIsStereo(newDetector);
117  break;
118  }
119 
120  const auto subdet = newDetector.subdetId();
121  const auto layer = tTopo.layer( newDetector );
122 
123  hasHit[subdet][layer] = true;
124  if(isPixel) hasPixel[subdet][layer] = isPixel;
125  else if(isStripMono) hasMono[subdet][layer] = isStripMono;
126  else if(isStripStereo) hasStereo[subdet][layer] = isStripStereo;
127  }
128  }
129  }
130 
131 
132  unsigned int nLayers = 0;
133  unsigned int nPixelLayers = 0;
134  unsigned int nStripMonoAndStereoLayers = 0;
135  for(unsigned int i=0; i<maxSubdet; ++i) {
136  for(unsigned int j=0; j<maxLayer; ++j) {
137  nLayers += hasHit[i][j];
138  nPixelLayers += hasPixel[i][j];
139  nStripMonoAndStereoLayers += (hasMono[i][j] && hasStereo[i][j]);
140  }
141  }
142 
143  valuesLayers[iTP] = nLayers;
144  valuesPixelLayers[iTP] = nPixelLayers;
145  valuesStripMonoAndStereoLayers[iTP] = nStripMonoAndStereoLayers;
146  }
147 
148  auto ret0 = std::make_unique<edm::ValueMap<unsigned int>>();
149  {
151  filler.insert(htps, valuesLayers.begin(), valuesLayers.end());
152  filler.fill();
153  }
154  auto ret1 = std::make_unique<edm::ValueMap<unsigned int>>();
155  {
157  filler.insert(htps, valuesPixelLayers.begin(), valuesPixelLayers.end());
158  filler.fill();
159  }
160  auto ret2 = std::make_unique<edm::ValueMap<unsigned int>>();
161  {
163  filler.insert(htps, valuesStripMonoAndStereoLayers.begin(), valuesStripMonoAndStereoLayers.end());
164  filler.fill();
165  }
166 
167  return std::make_tuple(std::move(ret0), std::move(ret1), std::move(ret2));
168 }
const std::vector< SimTrack > & g4Tracks() const
std::vector< TrackingParticle > TrackingParticleCollection
bool tobIsStereo(const DetId &id) const
int pdgId() const
PDG ID.
bool tobIsRPhi(const DetId &id) const
#define constexpr
bool tidIsStereo(const DetId &id) const
bool tecIsStereo(const DetId &id) const
bool tibIsRPhi(const DetId &id) const
simTrack
per collection params
EncodedEventId eventId() const
Definition: PSimHit.h:105
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:38
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
const T & get() const
Definition: EventSetup.h:59
bool tibIsStereo(const DetId &id) const
unsigned short processType() const
Definition: PSimHit.h:118
unsigned int layer(const DetId &id) const
EncodedEventId eventId() const
Signal source, crossing number.
int particleType() const
Definition: PSimHit.h:85
bool isPixel(HitType hitType)
Monte Carlo truth information used for tracking validation.
def move(src, dest)
Definition: eostools.py:510
unsigned int detUnitId() const
Definition: PSimHit.h:93

Member Data Documentation

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

Definition at line 60 of file TrackingParticleNumberOfLayers.h.

Referenced by calculate(), and TrackingParticleNumberOfLayers().