Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:24:08

0001 #ifndef PhysicsTools_RecoUtils_CheckHitPattern_H
0002 #define PhysicsTools_RecoUtils_CheckHitPattern_H
0003 
0004 /*
0005  * Determine if a track has hits in front of its assumed production point.
0006  * Also determine if it misses hits between its assumed production point and its innermost hit.
0007  *
0008  * FIXME: as it stands it is pretty inefficient for numerous reasons
0009  *        if used seriously it needs to be optimized and used properly... 
0010  */
0011 
0012 #include "DataFormats/TrackReco/interface/Track.h"
0013 #include "RecoVertex/VertexPrimitives/interface/VertexState.h"
0014 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
0015 #include "DataFormats/TrackReco/interface/HitPattern.h"
0016 #include <utility>
0017 #include <map>
0018 
0019 class DetId;
0020 
0021 class TrackerTopology;
0022 class TrackerGeometry;
0023 
0024 class CheckHitPattern {
0025 public:
0026   struct Result {
0027     // Number of hits track has in front of the vertex.
0028     unsigned int hitsInFrontOfVert;
0029     // Number of missing hits between the vertex position and the innermost valid hit on the track.
0030     unsigned int missHitsAfterVert;
0031   };
0032 
0033   // Check if hit pattern of this track is consistent with it being produced
0034   // at given vertex. See comments above for "Result" struct for details of returned information.
0035   Result operator()(const reco::Track& track, const VertexState& vert) const;
0036 
0037   // Print hit pattern on track
0038   static void print(const reco::Track& track);
0039 
0040   // Create map indicating r/z values of all layers/disks.
0041   void init(const TrackerTopology* tTopo, const TrackerGeometry& geom, const TransientTrackBuilder& builder);
0042 
0043   // Return a pair<uint32, uint32> consisting of the numbers used by HitPattern to
0044   // identify subdetector and layer number respectively.
0045   typedef std::pair<uint32_t, uint32_t> DetInfo;
0046   static DetInfo interpretDetId(DetId detId, const TrackerTopology* tTopo);
0047 
0048   // Return a bool indicating if a given subdetector is in the barrel.
0049   static bool barrel(uint32_t subDet);
0050 
0051   static void print(const reco::HitPattern::HitCategory category, const reco::HitPattern& hp);
0052 
0053 private:
0054   // Note if geometry info is already initialized.
0055   bool geomInitDone_ = false;
0056 
0057   // For a given subdetector & layer number, this stores the minimum and maximum
0058   // r (or z) values if it is barrel (or endcap) respectively.
0059   typedef std::map<DetInfo, std::pair<double, double> > RZrangeMap;
0060   RZrangeMap rangeRorZ_;
0061 
0062   // Makes TransientTracks needed for vertex fitting.
0063   const TransientTrackBuilder* trkTool_ = nullptr;
0064 };
0065 
0066 #endif