CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFCheckHitPattern.cc
Go to the documentation of this file.
2 
3 // To get Tracker Geometry
9 
10 // To convert detId to subdet/layer number.
19 
20 #include <map>
21 
22 using namespace reco;
23 using namespace std;
24 
26 
27  //
28  // Note min/max radius (z) of each barrel layer (endcap disk).
29  //
30 
31  geomInitDone_ = true;
32 
33  // Get Tracker geometry
34  const TrackingGeometry::DetContainer& dets = tkerGeomHandle_->dets();
35 
36  // Loop over all modules in the Tracker.
37  for (unsigned int i = 0; i < dets.size(); i++) {
38 
39  // Get subdet and layer of this module
40  DetInfo detInfo = this->interpretDetId(dets[i]->geographicalId());
41  uint32_t subDet = detInfo.first;
42 
43  // Note r (or z) of module if barrel (or endcap).
44  double r_or_z;
45  if (this->barrel(subDet)) {
46  r_or_z = dets[i]->position().perp();
47  } else {
48  r_or_z = fabs(dets[i]->position().z());
49  }
50 
51  // Recover min/max r/z value of this layer/disk found so far.
52  double minRZ = 999.;
53  double maxRZ = 0.;
54  if (rangeRorZ_.find(detInfo) != rangeRorZ_.end()) {
55  minRZ = rangeRorZ_[detInfo].first;
56  maxRZ = rangeRorZ_[detInfo].second;
57  }
58 
59  // Update with those of this module if it exceeds them.
60  if (minRZ > r_or_z) minRZ = r_or_z;
61  if (maxRZ < r_or_z) maxRZ = r_or_z;
62  rangeRorZ_[detInfo] = pair<double, double>(minRZ, maxRZ);
63  }
64 
65 #if 0
66  //def DEBUG_CHECKHITPATTERN
67  RZrangeMap::const_iterator d;
68  for (d = rangeRorZ_.begin(); d != rangeRorZ_.end(); d++) {
69  DetInfo detInfo = d->first;
70  pair<double, double> rangeRZ = d->second;
71  }
72 #endif
73 }
74 
76  // Convert detId to a pair<uint32, uint32> consisting of the numbers used by HitPattern
77  // to identify subdetector and layer number respectively.
78  if (detId.subdetId() == StripSubdetector::TIB) {
79  return DetInfo( detId.subdetId() , TIBDetId(detId).layer() );
80  } else if (detId.subdetId() == StripSubdetector::TOB) {
81  return DetInfo( detId.subdetId() , TOBDetId(detId).layer() );
82  } else if (detId.subdetId() == StripSubdetector::TID) {
83  return DetInfo( detId.subdetId() , TIDDetId(detId).wheel() );
84  } else if (detId.subdetId() == StripSubdetector::TEC) {
85  return DetInfo( detId.subdetId() , TECDetId(detId).wheel() );
86  } else if (detId.subdetId() == PixelSubdetector::PixelBarrel) {
87  return DetInfo( detId.subdetId() , PXBDetId(detId).layer() );
88  } else if (detId.subdetId() == PixelSubdetector::PixelEndcap) {
89  return DetInfo( detId.subdetId() , PXFDetId(detId).disk() );
90  } else {
91  throw cms::Exception("RecoParticleFlow", "Found DetId that is not in Tracker");
92  }
93 }
94 
95 bool PFCheckHitPattern::barrel(uint32_t subDet) {
96  // Determines if given sub-detector is in the barrel.
97  return (subDet == StripSubdetector::TIB || subDet == StripSubdetector::TOB ||
99 }
100 
101 
102 pair< PFCheckHitPattern::PFTrackHitInfo, PFCheckHitPattern::PFTrackHitInfo>
104  const TrackBaseRef track, const TransientVertex& vert)
105 {
106 
107  // PFCheck if hit pattern of this track is consistent with it being produced
108  // at given vertex. Pair.first gives number of hits on track in front of vertex.
109  // Pair.second gives number of missing hits between vertex and innermost hit
110  // on track.
111 
112  // Initialise geometry info if not yet done.
113  if (!geomInitDone_) this->init(tkerGeomHandle_);
114 
115  // Get hit patterns of this track
116  const reco::HitPattern& hp = track.get()->hitPattern();
117  const reco::HitPattern& ip = track.get()->trackerExpectedHitsInner();
118 
119  // Count number of valid hits on track definately in front of the vertex,
120  // taking into account finite depth of each layer.
121  unsigned int nHitBefore = 0;
122  unsigned int nHitAfter = 0;
123 
124  for (int i = 0; i < hp.numberOfHits(); i++) {
125  uint32_t hit = hp.getHitPattern(i);
126  if (hp.trackerHitFilter(hit) && hp.validHitFilter(hit)) {
127  uint32_t subDet = hp.getSubStructure(hit);
128  uint32_t layer = hp.getLayer(hit);
129  DetInfo detInfo(subDet, layer);
130  double maxRZ = rangeRorZ_[detInfo].second;
131 
132  if (this->barrel(subDet)) {
133  if (vert.position().perp() > maxRZ) nHitBefore++;
134  else nHitAfter++;
135  } else {
136  if (fabs(vert.position().z()) > maxRZ) nHitBefore++;
137  else nHitAfter++;
138  }
139  }
140  }
141 
142  // Count number of missing hits before the innermost hit on the track,
143  // taking into account finite depth of each layer.
144  unsigned int nMissHitAfter = 0;
145  unsigned int nMissHitBefore = 0;
146 
147  for (int i = 0; i < ip.numberOfHits(); i++) {
148  uint32_t hit = ip.getHitPattern(i);
149  if (ip.trackerHitFilter(hit)) {
150  uint32_t subDet = ip.getSubStructure(hit);
151  uint32_t layer = ip.getLayer(hit);
152  DetInfo detInfo(subDet, layer);
153  double minRZ = rangeRorZ_[detInfo].first;
154 
155  // cout << "subDet = " << subDet << " layer = " << layer << " minRZ = " << minRZ << endl;
156 
157  if (this->barrel(subDet)) {
158  if (vert.position().perp() < minRZ) nMissHitAfter++;
159  else nMissHitBefore++;
160  } else {
161  if (fabs(vert.position().z()) < minRZ) nMissHitAfter++;
162  else nMissHitBefore++;
163  }
164  }
165  }
166 
167 
168  PFTrackHitInfo trackToVertex(nHitBefore, nMissHitBefore);
169  PFTrackHitInfo trackFromVertex(nHitAfter, nMissHitAfter);
170 
171 
172  return pair< PFTrackHitInfo, PFTrackHitInfo>(trackToVertex, trackFromVertex);
173 }
174 
175 void PFCheckHitPattern::print(const TrackBaseRef track) const {
176  // Get hit patterns of this track
177  const reco::HitPattern& hp = track.get()->hitPattern();
178  const reco::HitPattern& ip = track.get()->trackerExpectedHitsInner();
179 
180  cout<<"=== Hits on Track ==="<<endl;
181  this->print(hp);
182  cout<<"=== Hits before track ==="<<endl;
183  this->print(ip);
184 }
185 
187  for (int i = 0; i < hp.numberOfHits(); i++) {
188  uint32_t hit = hp.getHitPattern(i);
189  if (hp.trackerHitFilter(hit)) {
190  uint32_t subdet = hp.getSubStructure(hit);
191  uint32_t layer = hp.getLayer(hit);
192  cout<<"hit "<<i<<" subdet="<<subdet<<" layer="<<layer<<" type "<<hp.getHitType(hit)<<endl;
193  }
194  }
195 }
int i
Definition: DBlmapReader.cc:9
T perp() const
Definition: PV3DBase.h:72
unsigned int layer() const
layer id
Definition: TOBDetId.h:39
static uint32_t getLayer(uint32_t pattern)
Definition: HitPattern.h:520
static bool trackerHitFilter(uint32_t pattern)
Definition: HitPattern.h:501
std::string print(const Track &, edm::Verbosity=edm::Concise)
Track print utility.
Definition: print.cc:8
int init
Definition: HydjetWrapper.h:62
static uint32_t getHitType(uint32_t pattern)
Definition: HitPattern.h:536
void print(const reco::TrackBaseRef track) const
Print hit pattern on track.
float float float z
unsigned int layer() const
layer id
Definition: PXBDetId.h:35
PFTrackHitFullInfo analyze(edm::ESHandle< TrackerGeometry >, const reco::TrackBaseRef track, const TransientVertex &vert)
std::pair< unsigned int, unsigned int > PFTrackHitInfo
static DetInfo interpretDetId(DetId detId)
susybsm::HSCParticleRefProd hp
Definition: classes.h:27
GlobalPoint position() const
T z() const
Definition: PV3DBase.h:64
int numberOfHits() const
Definition: HitPattern.cc:211
static uint32_t getSubStructure(uint32_t pattern)
Definition: HitPattern.h:514
const HitPattern & trackerExpectedHitsInner() const
Access the hit pattern counting (in the Tracker) the number of expected crossed layers before the fir...
Definition: TrackBase.h:223
static bool validHitFilter(uint32_t pattern)
Definition: HitPattern.h:564
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:221
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
unsigned int disk() const
disk id
Definition: PXFDetId.h:43
Definition: DetId.h:18
unsigned int wheel() const
wheel id
Definition: TECDetId.h:52
unsigned int layer() const
layer id
Definition: TIBDetId.h:41
std::pair< uint32_t, uint32_t > DetInfo
static int position[264][3]
Definition: ReadPGInfo.cc:509
void init(edm::ESHandle< TrackerGeometry >)
Create map indicating r/z values of all layers/disks.
tuple cout
Definition: gather_cfg.py:121
uint32_t getHitPattern(int position) const
Definition: HitPattern.cc:142
std::vector< GeomDet const * > DetContainer
value_type const * get() const
Definition: RefToBase.h:212
static bool barrel(uint32_t subDet)
Return a bool indicating if a given subdetector is in the barrel.
unsigned int wheel() const
wheel id
Definition: TIDDetId.h:50