CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 14 of file TrackingParticleNumberOfLayers.cc.

References edm::Event::getByToken(), and trackIdToHitPtr_.

14  {
15 
16  // A multimap linking SimTrack::trackId() to a pointer to PSimHit
17  // Similar to TrackingTruthAccumulator
18  for(const auto& simHitToken: simHitTokens) {
20  iEvent.getByToken(simHitToken, hsimhits);
21  trackIdToHitPtr_.reserve(trackIdToHitPtr_.size()+hsimhits->size());
22  for(const auto& simHit: *hsimhits) {
23  trackIdToHitPtr_.emplace_back(simHit.trackId(), &simHit);
24  }
25  }
26  std::stable_sort(trackIdToHitPtr_.begin(), trackIdToHitPtr_.end(), trackIdHitPairLess);
27 }
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
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 32 of file TrackingParticleNumberOfLayers.cc.

References constexpr, PSimHit::detUnitId(), TrackingParticle::eventId(), PSimHit::eventId(), edm::helper::Filler< Map >::fill(), TrackingParticle::g4Tracks(), edm::EventSetup::get(), i, edm::helper::Filler< Map >::insert(), fastTrackerRecHitType::isPixel(), j, TrackerTopology::layer(), eostools::move(), nPixelLayers, nStripMonoAndStereoLayers, PSimHit::particleType(), HLT_25ns10e33_v2_cff::particleType, TrackingParticle::pdgId(), SingleMuPt40Fwdv2_cfi_GEN_SIM::pdgId, PixelSubdetector::PixelBarrel, PixelSubdetector::PixelEndcap, PSimHit::processType(), 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 MultiTrackValidator::analyze(), TrackerSeedValidator::analyze(), and TrackingParticleNumberOfLayersProducer::produce().

32  {
33  //Retrieve tracker topology from geometry
35  iSetup.get<TrackerTopologyRcd>().get(tTopoHandle);
36  const TrackerTopology& tTopo = *tTopoHandle;
37 
38  const TrackingParticleCollection& tps = *htps;
39  std::vector<unsigned int> valuesLayers(tps.size(), 0);
40  std::vector<unsigned int> valuesPixelLayers(tps.size(), 0);
41  std::vector<unsigned int> valuesStripMonoAndStereoLayers(tps.size(), 0);
42  for(size_t iTP=0; iTP<tps.size(); ++iTP) {
43  const TrackingParticle& tp = tps[iTP];
44  const auto pdgId = tp.pdgId();
45 
46  // I would prefer a better way...
47  constexpr unsigned int maxSubdet = static_cast<unsigned>(StripSubdetector::TEC)+1;
48  constexpr unsigned int maxLayer = 0xF+1; // as in HitPattern.h
49  bool hasHit[maxSubdet][maxLayer];
50  bool hasPixel[maxSubdet][maxLayer];
51  bool hasMono[maxSubdet][maxLayer];
52  bool hasStereo[maxSubdet][maxLayer];
53  memset(hasHit, 0, sizeof(hasHit));
54  memset(hasPixel, 0, sizeof(hasPixel));
55  memset(hasMono, 0, sizeof(hasMono));
56  memset(hasStereo, 0, sizeof(hasStereo));
57 
58  for(const SimTrack& simTrack: tp.g4Tracks()) {
59  // Logic is from TrackingTruthAccumulator
60  auto range = std::equal_range(trackIdToHitPtr_.begin(), trackIdToHitPtr_.end(), std::pair<unsigned int, const PSimHit *>(simTrack.trackId(), nullptr), trackIdHitPairLess);
61  if(range.first == range.second) continue;
62 
63  auto iHitPtr = range.first;
64  int processType = iHitPtr->second->processType();
65  int particleType = iHitPtr->second->particleType();
66 
67  for(; iHitPtr != range.second; ++iHitPtr) {
68  const PSimHit& simHit = *(iHitPtr->second);
69  if(simHit.eventId() != tp.eventId())
70  continue;
71  DetId newDetector = DetId( simHit.detUnitId() );
72 
73  // Check for delta and interaction products discards
74  if( processType == simHit.processType() && particleType == simHit.particleType() && pdgId == simHit.particleType() ) {
75  // The logic of this piece follows HitPattern
76  bool isPixel = false;
77  bool isStripStereo = false;
78  bool isStripMono = false;
79 
80  switch(newDetector.subdetId()) {
83  isPixel = true;
84  break;
86  isStripMono = tTopo.tibIsRPhi(newDetector);
87  isStripStereo = tTopo.tibIsStereo(newDetector);
88  break;
90  isStripMono = tTopo.tidIsRPhi(newDetector);
91  isStripStereo = tTopo.tidIsStereo(newDetector);
92  break;
94  isStripMono = tTopo.tobIsRPhi(newDetector);
95  isStripStereo = tTopo.tobIsStereo(newDetector);
96  break;
98  isStripMono = tTopo.tecIsRPhi(newDetector);
99  isStripStereo = tTopo.tecIsStereo(newDetector);
100  break;
101  }
102 
103  const auto subdet = newDetector.subdetId();
104  const auto layer = tTopo.layer( newDetector );
105 
106  hasHit[subdet][layer] = true;
107  if(isPixel) hasPixel[subdet][layer] = isPixel;
108  else if(isStripMono) hasMono[subdet][layer] = isStripMono;
109  else if(isStripStereo) hasStereo[subdet][layer] = isStripStereo;
110  }
111  }
112  }
113 
114 
115  unsigned int nLayers = 0;
116  unsigned int nPixelLayers = 0;
117  unsigned int nStripMonoAndStereoLayers = 0;
118  for(unsigned int i=0; i<maxSubdet; ++i) {
119  for(unsigned int j=0; j<maxLayer; ++j) {
120  nLayers += hasHit[i][j];
121  nPixelLayers += hasPixel[i][j];
122  nStripMonoAndStereoLayers += (hasMono[i][j] && hasStereo[i][j]);
123  }
124  }
125 
126  valuesLayers[iTP] = nLayers;
127  valuesPixelLayers[iTP] = nPixelLayers;
128  valuesStripMonoAndStereoLayers[iTP] = nStripMonoAndStereoLayers;
129  }
130 
131  auto ret0 = std::make_unique<edm::ValueMap<unsigned int>>();
132  {
134  filler.insert(htps, valuesLayers.begin(), valuesLayers.end());
135  filler.fill();
136  }
137  auto ret1 = std::make_unique<edm::ValueMap<unsigned int>>();
138  {
140  filler.insert(htps, valuesPixelLayers.begin(), valuesPixelLayers.end());
141  filler.fill();
142  }
143  auto ret2 = std::make_unique<edm::ValueMap<unsigned int>>();
144  {
146  filler.insert(htps, valuesStripMonoAndStereoLayers.begin(), valuesStripMonoAndStereoLayers.end());
147  filler.fill();
148  }
149 
150  return std::make_tuple(std::move(ret0), std::move(ret1), std::move(ret2));
151 }
int i
Definition: DBlmapReader.cc:9
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
def move
Definition: eostools.py:510
int j
Definition: DBlmapReader.cc:9
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:37
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:56
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.
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().