CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CheckHitPattern.h
Go to the documentation of this file.
1 #ifndef PhysicsTools_RecoUtils_CheckHitPattern_H
2 #define PhysicsTools_RecoUtils_CheckHitPattern_H
3 
4 /*
5  * Determine if a track has hits in front of its assumed production point.
6  * Also determine if it misses hits between its assumed production point and its innermost hit.
7  */
8 
9 // standard EDAnalyser include files
10 #include <memory>
17 
22 
23 #include <utility>
24 #include <map>
25 
26 class DetId;
27 
29 
30 public:
31 
32  struct Result {
33  // Number of hits track has in front of the vertex.
34  unsigned int hitsInFrontOfVert;
35  // Number of missing hits between the vertex position and the innermost valid hit on the track.
36  unsigned int missHitsAfterVert;
37  };
38 
40 
42 
43  // Check if hit pattern of this track is consistent with it being produced
44  // at given vertex. See comments above for "Result" struct for details of returned information.
45  // N.B. If FixHitPattern = true, then Result.missHitsAfterVert will be calculated after rederiving
46  // the missing hit pattern. This rederivation is sometimes a good idea, since otherwise the
47  // number of missing hits can be substantially underestimated. See comments in FixTrackHitPattern.h
48  // for details.
49  Result analyze(const edm::EventSetup& iSetup,
50  const reco::Track& track, const VertexState& vert, bool fixHitPattern=true);
51 
52  // Print hit pattern on track
53  void print(const reco::Track& track) const;
54 
55 private:
56  // Create map indicating r/z values of all layers/disks.
57  void init (const edm::EventSetup& iSetup);
58 
59  // Return a pair<uint32, uint32> consisting of the numbers used by HitPattern to
60  // identify subdetector and layer number respectively.
61  typedef std::pair<uint32_t, uint32_t> DetInfo;
62  static DetInfo interpretDetId(DetId detId);
63 
64  // Return a bool indicating if a given subdetector is in the barrel.
65  static bool barrel(uint32_t subDet);
66 
67  void print(const reco::HitPattern& hp) const;
68 
69 private:
70  // Note if geometry info is already initialized.
72 
73  // For a given subdetector & layer number, this stores the minimum and maximum
74  // r (or z) values if it is barrel (or endcap) respectively.
75  typedef std::map< DetInfo, std::pair< double, double> > RZrangeMap;
77 
78  // Makes TransientTracks needed for vertex fitting.
80 };
81 
82 #endif
std::map< DetInfo, std::pair< double, double > > RZrangeMap
Result analyze(const edm::EventSetup &iSetup, const reco::Track &track, const VertexState &vert, bool fixHitPattern=true)
std::pair< uint32_t, uint32_t > DetInfo
static DetInfo interpretDetId(DetId detId)
unsigned int missHitsAfterVert
unsigned int hitsInFrontOfVert
Definition: DetId.h:20
void init(const edm::EventSetup &iSetup)
static RZrangeMap rangeRorZ_
edm::ESHandle< TransientTrackBuilder > trkTool_
static bool barrel(uint32_t subDet)
void print(const reco::Track &track) const