CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CheckHitPattern.cc
Go to the documentation of this file.
3 
4 // To get Tracker Geometry
8 
9 // To convert detId to subdet/layer number.
12 //#include "DataFormats/SiStripDetId/interface/SiStripDetId.h"
19 
22 
24 
25 // For a given subdetector & layer number, this static map stores the minimum and maximum
26 // r (or z) values if it is barrel (or endcap) respectively.
28 
30 
31  //
32  // Note min/max radius (z) of each barrel layer (endcap disk).
33  //
34 
35  geomInitDone_ = true;
36 
37  // Get Tracker geometry
38  edm::ESHandle<TrackerGeometry> trackerGeometry;
39  iSetup.get<TrackerDigiGeometryRecord>().get(trackerGeometry);
40  const TrackingGeometry::DetContainer& dets = trackerGeometry->dets();
41 
42  // Loop over all modules in the Tracker.
43  for (unsigned int i = 0; i < dets.size(); i++) {
44 
45  // Get subdet and layer of this module
46  DetInfo detInfo = this->interpretDetId(dets[i]->geographicalId());
47  uint32_t subDet = detInfo.first;
48 
49  // Note r (or z) of module if barrel (or endcap).
50  double r_or_z;
51  if (this->barrel(subDet)) {
52  r_or_z = dets[i]->position().perp();
53  } else {
54  r_or_z = fabs(dets[i]->position().z());
55  }
56 
57  // Recover min/max r/z value of this layer/disk found so far.
58  double minRZ = 999.;
59  double maxRZ = 0.;
60  if (rangeRorZ_.find(detInfo) != rangeRorZ_.end()) {
61  minRZ = rangeRorZ_[detInfo].first;
62  maxRZ = rangeRorZ_[detInfo].second;
63  }
64 
65  // Update with those of this module if it exceeds them.
66  if (minRZ > r_or_z) minRZ = r_or_z;
67  if (maxRZ < r_or_z) maxRZ = r_or_z;
68  rangeRorZ_[detInfo] = std::pair<double, double>(minRZ, maxRZ);
69  }
70 
71 #ifdef DEBUG_CHECKHITPATTERN
72  RZrangeMap::const_iterator d;
73  for (d = rangeRorZ_.begin(); d != rangeRorZ_.end(); d++) {
74  DetInfo detInfo = d->first;
75  std::pair<double, double> rangeRZ = d->second;
76  std::std::cout<<"CHECKHITPATTERN: Tracker subdetector type="<<detInfo.first<<" layer="<<detInfo.second
77  <<" has min r (or z) ="<<rangeRZ.first<<" and max r (or z) = "<<rangeRZ.second<<std::std::endl;
78  }
79 #endif
80 }
81 
83  // Convert detId to a pair<uint32, uint32> consisting of the numbers used by HitPattern
84  // to identify subdetector and layer number respectively.
85  if (detId.subdetId() == StripSubdetector::TIB) {
86  return DetInfo( detId.subdetId() , TIBDetId(detId).layer() );
87  } else if (detId.subdetId() == StripSubdetector::TOB) {
88  return DetInfo( detId.subdetId() , TOBDetId(detId).layer() );
89  } else if (detId.subdetId() == StripSubdetector::TID) {
90  return DetInfo( detId.subdetId() , TIDDetId(detId).wheel() );
91  } else if (detId.subdetId() == StripSubdetector::TEC) {
92  return DetInfo( detId.subdetId() , TECDetId(detId).wheel() );
93  } else if (detId.subdetId() == PixelSubdetector::PixelBarrel) {
94  return DetInfo( detId.subdetId() , PXBDetId(detId).layer() );
95  } else if (detId.subdetId() == PixelSubdetector::PixelEndcap) {
96  return DetInfo( detId.subdetId() , PXFDetId(detId).disk() );
97  } else {
98  throw cms::Exception("NotFound","Found DetId that is not in Tracker");
99  }
100 }
101 
102 bool CheckHitPattern::barrel(uint32_t subDet) {
103  // Determines if given sub-detector is in the barrel.
104  return (subDet == StripSubdetector::TIB || subDet == StripSubdetector::TOB ||
105  subDet == PixelSubdetector::PixelBarrel);
106 }
107 
108 
110  const reco::Track& track, const VertexState& vert, bool fixHitPattern)
111 {
112  // Check if hit pattern of this track is consistent with it being produced
113  // at given vertex.
114 
115  // Initialise geometry info if not yet done.
116  if (!geomInitDone_) this->init(iSetup);
117 
118  // Optionally set vertex position to zero for debugging.
119  // VertexState vertDebug( GlobalPoint(0.,0.,0.) , GlobalError(1e-8, 0., 1e-8, 0., 0., 1e-8) );
120 
121  // Evaluate track parameters at vertex.
122  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",trkTool_); // Needed for vertex fits
123  reco::TransientTrack t_trk = trkTool_->build(track);
124  GlobalVector p3_trk = t_trk.trajectoryStateClosestToPoint(vert.position()).momentum();
125  bool trkGoesInsideOut = fabs(reco::deltaPhi<const GlobalVector, const GlobalPoint>(p3_trk, vert.position())) < 0.5*M_PI;
126 
127  LogDebug("CHP")<<"TRACK: in-->out ? "<<trkGoesInsideOut<<" dxy="<<track.dxy()<<" sz="<<track.dz()<<" eta="<<track.eta()<<" barrel hits="<<track.hitPattern().numberOfValidPixelHits()<<"/"<<track.hitPattern().numberOfValidStripTIBHits()<<"/"<<track.hitPattern().numberOfValidStripTOBHits();
128  LogDebug("CHP")<<"VERT: r="<<vert.position().perp()<<" z="<<vert.position().z();
129  // 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";
130  // Get hit patterns of this track
131  const reco::HitPattern& hp = track.hitPattern();
133 
134  // Optionally fix inner hit pattern (needed if uncertainty on track trajectory is large).
135  if (fixHitPattern) {
136  static FixTrackHitPattern fixTrackHitPattern;
137  FixTrackHitPattern::Result fixedHP = fixTrackHitPattern.analyze(iSetup, track);
138  ip = fixedHP.innerHitPattern;
139  }
140 
141  // Count number of valid hits on track definately in front of the vertex,
142  // taking into account finite depth of each layer.
143  unsigned int nHitBefore = 0;
144  for (int i = 0; i < hp.numberOfHits(); i++) {
145  uint32_t hit = hp.getHitPattern(i);
146  if (hp.trackerHitFilter(hit) && hp.validHitFilter(hit)) {
147  uint32_t subDet = hp.getSubStructure(hit);
148  uint32_t layer = hp.getLayer(hit);
149  DetInfo detInfo(subDet, layer);
150  double maxRZ = rangeRorZ_[detInfo].second;
151 
152  if (this->barrel(subDet)) {
153  // Be careful. If the track starts by going outside-->in, it is allowed to have hits before the vertex !
154  if (vert.position().perp() > maxRZ && trkGoesInsideOut) nHitBefore++;
155  } else {
156  if (fabs(vert.position().z()) > maxRZ) nHitBefore++;
157  }
158  }
159  }
160 
161  // Count number of missing hits before the innermost hit on the track,
162  // taking into account finite depth of each layer.
163  unsigned int nMissHitAfter = 0;
164  for (int i = 0; i < ip.numberOfHits(); i++) {
165  uint32_t hit = ip.getHitPattern(i);
166  // if (ip.trackerHitFilter(hit)) {
167  if (ip.trackerHitFilter(hit) && ip.type_1_HitFilter(hit)) {
168  uint32_t subDet = ip.getSubStructure(hit);
169  uint32_t layer = ip.getLayer(hit);
170  DetInfo detInfo(subDet, layer);
171  double minRZ = rangeRorZ_[detInfo].first;
172 
173  if (this->barrel(subDet)) {
174  // Be careful. If the track starts by going outside-->in, then it misses hits
175  // in all layers it crosses before its innermost valid hit.
176  if (vert.position().perp() < minRZ || ! trkGoesInsideOut) nMissHitAfter++;
177  } else {
178  if (fabs(vert.position().z()) < minRZ) nMissHitAfter++;
179  }
180  }
181  }
182 
183  Result result;
184  result.hitsInFrontOfVert = nHitBefore;
185  result.missHitsAfterVert = nMissHitAfter;
186  return result;
187 }
188 
189 void CheckHitPattern::print(const reco::Track& track) const {
190  // Get hit patterns of this track
191  const reco::HitPattern& hp = track.hitPattern();
192  const reco::HitPattern& ip = track.trackerExpectedHitsInner();
193 
194  std::cout<<"=== Hits on Track ==="<<std::endl;
195  this->print(hp);
196  std::cout<<"=== Hits before track ==="<<std::endl;
197  this->print(ip);
198 }
199 
201  for (int i = 0; i < hp.numberOfHits(); i++) {
202  uint32_t hit = hp.getHitPattern(i);
203  if (hp.trackerHitFilter(hit)) {
204  uint32_t subdet = hp.getSubStructure(hit);
205  uint32_t layer = hp.getLayer(hit);
206  std::cout<<"hit "<<i<<" subdet="<<subdet<<" layer="<<layer<<" type "<<hp.getHitType(hit)<<std::endl;
207  }
208  }
209 }
#define LogDebug(id)
int i
Definition: DBlmapReader.cc:9
std::map< DetInfo, std::pair< double, double > > RZrangeMap
T perp() const
Definition: PV3DBase.h:71
unsigned int layer() const
layer id
Definition: TOBDetId.h:39
Result analyze(const edm::EventSetup &iSetup, const reco::Track &track, const VertexState &vert, bool fixHitPattern=true)
static uint32_t getLayer(uint32_t pattern)
Definition: HitPattern.h:485
std::pair< uint32_t, uint32_t > DetInfo
static DetInfo interpretDetId(DetId detId)
static bool trackerHitFilter(uint32_t pattern)
Definition: HitPattern.h:466
static uint32_t getHitType(uint32_t pattern)
Definition: HitPattern.h:501
std::vector< GeomDet * > DetContainer
GlobalPoint position() const
Definition: VertexState.h:29
int numberOfValidStripTOBHits() const
Definition: HitPattern.h:590
double double double z
static int position[TOTALCHAMBERS][3]
Definition: ReadPGInfo.cc:509
unsigned int layer() const
layer id
Definition: PXBDetId.h:35
unsigned int missHitsAfterVert
unsigned int hitsInFrontOfVert
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:141
Result analyze(const edm::EventSetup &iSetup, const reco::Track &track)
static bool type_1_HitFilter(uint32_t pattern)
Definition: HitPattern.h:534
T z() const
Definition: PV3DBase.h:63
tuple result
Definition: query.py:137
int numberOfHits() const
Definition: HitPattern.cc:213
static uint32_t getSubStructure(uint32_t pattern)
Definition: HitPattern.h:479
const HitPattern & trackerExpectedHitsInner() const
Access the hit pattern counting (in the Tracker) the number of expected crossed layers before the fir...
Definition: TrackBase.h:225
static bool validHitFilter(uint32_t pattern)
Definition: HitPattern.h:529
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:223
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:39
unsigned int disk() const
disk id
Definition: PXFDetId.h:43
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:127
Definition: DetId.h:20
void init(const edm::EventSetup &iSetup)
reco::HitPattern innerHitPattern
#define M_PI
Definition: BFit3D.cc:3
TrajectoryStateClosestToPoint trajectoryStateClosestToPoint(const GlobalPoint &point) const
const T & get() const
Definition: EventSetup.h:55
int numberOfValidStripTIBHits() const
Definition: HitPattern.h:582
unsigned int wheel() const
wheel id
Definition: TECDetId.h:52
static RZrangeMap rangeRorZ_
unsigned int layer() const
layer id
Definition: TIBDetId.h:41
int numberOfValidPixelHits() const
Definition: HitPattern.h:566
edm::ESHandle< TransientTrackBuilder > trkTool_
static bool barrel(uint32_t subDet)
tuple cout
Definition: gather_cfg.py:121
uint32_t getHitPattern(int position) const
Definition: HitPattern.cc:144
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:121
void print(const reco::Track &track) const
unsigned int wheel() const
wheel id
Definition: TIDDetId.h:50