CMS 3D CMS Logo

CheckHitPattern.cc
Go to the documentation of this file.
2 
3 // To get Tracker Geometry
6 
7 // To convert detId to subdet/layer number.
10 //#include "DataFormats/SiStripDetId/interface/SiStripDetId.h"
13 
16 
19 
21  const TrackerGeometry& geom,
22  const TransientTrackBuilder& builder) {
23  trkTool_ = &builder;
24  //
25  // Note min/max radius (z) of each barrel layer (endcap disk).
26  //
27 
28  geomInitDone_ = true;
29 
30  const TrackingGeometry::DetContainer& dets = geom.dets();
31 
32  // Loop over all modules in the Tracker.
33  for (unsigned int i = 0; i < dets.size(); i++) {
34  // Get subdet and layer of this module
35  DetInfo detInfo = this->interpretDetId(dets[i]->geographicalId(), tTopo);
36  uint32_t subDet = detInfo.first;
37 
38  // Note r (or z) of module if barrel (or endcap).
39  double r_or_z;
40  if (this->barrel(subDet)) {
41  r_or_z = dets[i]->position().perp();
42  } else {
43  r_or_z = fabs(dets[i]->position().z());
44  }
45 
46  // Recover min/max r/z value of this layer/disk found so far.
47  double minRZ = 999.;
48  double maxRZ = 0.;
49  if (rangeRorZ_.find(detInfo) != rangeRorZ_.end()) {
50  minRZ = rangeRorZ_[detInfo].first;
51  maxRZ = rangeRorZ_[detInfo].second;
52  }
53 
54  // Update with those of this module if it exceeds them.
55  if (minRZ > r_or_z)
56  minRZ = r_or_z;
57  if (maxRZ < r_or_z)
58  maxRZ = r_or_z;
59  rangeRorZ_[detInfo] = std::pair<double, double>(minRZ, maxRZ);
60  }
61 
62 #ifdef DEBUG_CHECKHITPATTERN
63  RZrangeMap::const_iterator d;
64  for (d = rangeRorZ_.begin(); d != rangeRorZ_.end(); d++) {
65  DetInfo detInfo = d->first;
66  std::pair<double, double> rangeRZ = d->second;
67  std::std::cout << "CHECKHITPATTERN: Tracker subdetector type=" << detInfo.first << " layer=" << detInfo.second
68  << " has min r (or z) =" << rangeRZ.first << " and max r (or z) = " << rangeRZ.second
69  << std::std::endl;
70  }
71 #endif
72 }
73 
75  // Convert detId to a pair<uint32, uint32> consisting of the numbers used by HitPattern
76  // to identify subdetector and layer number respectively.
77  return DetInfo(detId.subdetId(), tTopo->layer(detId));
78 }
79 
80 bool CheckHitPattern::barrel(uint32_t subDet) {
81  // Determines if given sub-detector is in the barrel.
82  return (subDet == StripSubdetector::TIB || subDet == StripSubdetector::TOB ||
84 }
85 
87  // Check if hit pattern of this track is consistent with it being produced
88  // at given vertex.
89 
90  // Initialise geometry info if not yet done.
91  if (!geomInitDone_)
92  throw cms::Exception("CheckHitPattern::operator() called before CheckHitPattern::init");
93 
94  // Optionally set vertex position to zero for debugging.
95  // VertexState vertDebug( GlobalPoint(0.,0.,0.) , GlobalError(1e-8, 0., 1e-8, 0., 0., 1e-8) );
96 
97  // Evaluate track parameters at vertex.
99  GlobalVector p3_trk = t_trk.trajectoryStateClosestToPoint(vert.position()).momentum();
100  bool trkGoesInsideOut =
101  fabs(reco::deltaPhi<const GlobalVector, const GlobalPoint>(p3_trk, vert.position())) < 0.5 * M_PI;
102 
103  LogDebug("CHP") << "TRACK: in-->out ? " << trkGoesInsideOut << " dxy=" << track.dxy() << " sz=" << track.dz()
104  << " eta=" << track.eta() << " barrel hits=" << track.hitPattern().numberOfValidPixelHits() << "/"
105  << track.hitPattern().numberOfValidStripTIBHits() << "/"
106  << track.hitPattern().numberOfValidStripTOBHits();
107  LogDebug("CHP") << "VERT: r=" << vert.position().perp() << " z=" << vert.position().z();
108  // if (vert.position().perp() < 3.5 && fabs(vert.position().z()) < 10. && fabs(track.eta()) < 1 && fabs(track.dxy()) < 2 && fabs(track.dz()) < 2 && track.hitPattern().numberOfValidPixelHits() == 0 && track.hitPattern().numberOfValidStripTIBHits() == 0) LogDebug("CHP")<<"LOOKATTHISTRACK";
109  // Get hit patterns of this track
110  const reco::HitPattern& hp = track.hitPattern();
111 
112  // Count number of valid hits on track definately in front of the vertex,
113  // taking into account finite depth of each layer.
114  unsigned int nHitBefore = 0;
115  for (int i = 0; i < hp.numberOfAllHits(reco::HitPattern::TRACK_HITS); i++) {
116  uint32_t hit = hp.getHitPattern(reco::HitPattern::TRACK_HITS, i);
118  uint32_t subDet = reco::HitPattern::getSubStructure(hit);
120  DetInfo detInfo(subDet, layer);
121  auto maxRZ = (*rangeRorZ_.find(detInfo)).second.second;
122 
123  if (this->barrel(subDet)) {
124  // Be careful. If the track starts by going outside-->in, it is allowed to have hits before the vertex !
125  if (vert.position().perp() > maxRZ && trkGoesInsideOut)
126  nHitBefore++;
127  } else {
128  if (fabs(vert.position().z()) > maxRZ)
129  nHitBefore++;
130  }
131  }
132  }
133 
134  // Count number of missing hits before the innermost hit on the track,
135  // taking into account finite depth of each layer.
136  unsigned int nMissHitAfter = 0;
137  for (int i = 0; i < hp.numberOfAllHits(reco::HitPattern::MISSING_INNER_HITS); i++) {
138  uint32_t hit = hp.getHitPattern(reco::HitPattern::MISSING_INNER_HITS, i);
139  // if (hp.trackerHitFilter(hit)) {
141  uint32_t subDet = reco::HitPattern::getSubStructure(hit);
143  DetInfo detInfo(subDet, layer);
144  auto minRZ = (*rangeRorZ_.find(detInfo)).second.first;
145 
146  if (this->barrel(subDet)) {
147  // Be careful. If the track starts by going outside-->in, then it misses hits
148  // in all layers it crosses before its innermost valid hit.
149  if (vert.position().perp() < minRZ || !trkGoesInsideOut)
150  nMissHitAfter++;
151  } else {
152  if (fabs(vert.position().z()) < minRZ)
153  nMissHitAfter++;
154  }
155  }
156  }
157 
158  Result result;
159  result.hitsInFrontOfVert = nHitBefore;
160  result.missHitsAfterVert = nMissHitAfter;
161  return result;
162 }
163 
165  // Get hit patterns of this track
166  const reco::HitPattern& hp = track.hitPattern();
167  std::cout << "=== Hits on Track ===" << std::endl;
169  std::cout << "=== Hits before track ===" << std::endl;
171 }
172 
174  for (int i = 0; i < hp.numberOfAllHits(category); i++) {
175  uint32_t hit = hp.getHitPattern(category, i);
177  uint32_t subdet = reco::HitPattern::getSubStructure(hit);
179  std::cout << "hit " << i << " subdet=" << subdet << " layer=" << layer << " type " << hp.getHitType(hit)
180  << std::endl;
181  }
182  }
183 }
static uint32_t getLayer(uint16_t pattern)
Definition: HitPattern.h:721
T perp() const
Definition: PV3DBase.h:69
const TransientTrackBuilder * trkTool_
void init(const TrackerTopology *tTopo, const TrackerGeometry &geom, const TransientTrackBuilder &builder)
std::pair< uint32_t, uint32_t > DetInfo
T z() const
Definition: PV3DBase.h:61
static bool missingHitFilter(uint16_t pattern)
Definition: HitPattern.h:798
unsigned int layer(const DetId &id) const
static bool validHitFilter(uint16_t pattern)
Definition: HitPattern.h:796
reco::TransientTrack build(const reco::Track *p) const
constexpr std::array< uint8_t, layerIndexSize< TrackerTraits > > layer
U second(std::pair< T, U > const &p)
std::vector< const GeomDet * > DetContainer
TrajectoryStateClosestToPoint trajectoryStateClosestToPoint(const GlobalPoint &point) 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
Result operator()(const reco::Track &track, const VertexState &vert) const
d
Definition: ztail.py:151
static uint32_t getSubStructure(uint16_t pattern)
Definition: HitPattern.h:713
#define M_PI
Definition: DetId.h:17
static bool trackerHitFilter(uint16_t pattern)
Definition: HitPattern.h:679
static constexpr auto TIB
static DetInfo interpretDetId(DetId detId, const TrackerTopology *tTopo)
static void print(const reco::Track &track)
static int position[264][3]
Definition: ReadPGInfo.cc:289
static bool barrel(uint32_t subDet)
RZrangeMap rangeRorZ_
#define LogDebug(id)
GlobalPoint position() const
Definition: VertexState.h:62