Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2021-02-14 12:50:02

0001 #ifndef RecoAlgos_ClusterStorer_h
0002 #define RecoAlgos_ClusterStorer_h
0003 /** \class ClusterStorer
0004  *
0005  * Helper to store clones of SiStrip- and SiPixelClusters
0006  * of selected RecHits
0007  * 
0008  * \author Gero Flucke, DESY
0009  *         based on code originally part of TrackSelector.h by Giovanni Petrucciani,
0010  *         but extended to deal with both SiStripRecHit1D and SiStripRecHit2D
0011  * 
0012  * \version $Revision: 1.1 $
0013  *
0014  * $Id: TrackSelector.h,v 1.1 2009/03/04 13:11:28 llista Exp $
0015  *
0016  */
0017 #include "DataFormats/TrackingRecHit/interface/TrackingRecHitFwd.h"
0018 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
0019 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2D.h"
0020 #include "DataFormats/TrackerRecHit2D/interface/Phase2TrackerRecHit1D.h"
0021 #include "DataFormats/SiStripCluster/interface/SiStripCluster.h"
0022 #include "DataFormats/SiPixelCluster/interface/SiPixelCluster.h"
0023 #include "DataFormats/Common/interface/DetSetVectorNew.h"
0024 
0025 namespace helper {
0026 
0027   class ClusterStorer {
0028   public:
0029     ClusterStorer() {}
0030     /// add cluster of newHit to list (throws if hit is of unknown type)
0031     void addCluster(TrackingRecHitCollection &hits, size_t index);
0032     /// clear records
0033     void clear();
0034     //------------------------------------------------------------------
0035     //!  Processes all the clusters of the tracks
0036     //!  (after the tracks have been dealt with),
0037     //!  need Refs to products (i.e. full collections) in the event.
0038     //------------------------------------------------------------------
0039     void processAllClusters(edmNew::DetSetVector<SiPixelCluster> &pixelDsvToFill,
0040                             edm::RefProd<edmNew::DetSetVector<SiPixelCluster> > refPixelClusters,
0041                             edmNew::DetSetVector<SiStripCluster> &stripDsvToFill,
0042                             edm::RefProd<edmNew::DetSetVector<SiStripCluster> > refStripClusters);
0043 
0044   private:
0045     /// A struct for clusters associated to hits
0046     template <typename ClusterRefType>
0047     class ClusterHitRecord {
0048     public:
0049       /// Create a record for a hit with a given index in the TrackingRecHitCollection.
0050       /// 'RecHitType' must have a method 'cluster()' that returns a 'ClusterRefType'.
0051       template <typename RecHitType>
0052       ClusterHitRecord(const RecHitType &hit, TrackingRecHitCollection &hits, size_t idx)
0053           : detid_(hit.geographicalId().rawId()), hits_(&hits), index_(idx), ref_(hit.cluster()) {}
0054       /// returns the detid
0055       uint32_t detid() const { return detid_; }
0056       /// this method is to be able to compare and see if two refs are the same
0057       const ClusterRefType &clusterRef() const { return ref_; }
0058       /// this one is to sort by detid and then by index of the rechit
0059       bool operator<(const ClusterHitRecord<ClusterRefType> &other) const {
0060         return (detid_ != other.detid_) ? detid_ < other.detid_ : ref_ < other.ref_;
0061       }
0062       /// Set the reference of the hit of this record to 'newRef',
0063       /// will not modify the ref stored in this object.
0064       template <typename RecHitType>
0065       void rekey(const ClusterRefType &newRef) const;
0066 
0067     private:
0068       ClusterHitRecord() {}  /// private => unusable
0069       uint32_t detid_;
0070       TrackingRecHitCollection *hits_;
0071       size_t index_;
0072       ClusterRefType ref_;
0073     };
0074 
0075     typedef ClusterHitRecord<SiPixelRecHit::ClusterRef> PixelClusterHitRecord;
0076     /// Assuming that the ClusterRef is the same for all SiStripRecHit*:
0077     typedef ClusterHitRecord<SiStripRecHit2D::ClusterRef> StripClusterHitRecord;
0078     //FIXME:: this is just temporary solution for phase2,
0079     //probably is good to add a Phase2ClusterStorer?
0080     typedef ClusterHitRecord<Phase2TrackerRecHit1D::CluRef> Phase2OTClusterHitRecord;
0081 
0082     //------------------------------------------------------------------
0083     //!  Processes all the clusters of a specific type
0084     //!  (after the tracks have been dealt with)
0085     //------------------------------------------------------------------
0086     template <typename HitType, typename ClusterType>
0087     void processClusters(std::vector<ClusterHitRecord<typename HitType::ClusterRef> > &clusterRecords,
0088                          edmNew::DetSetVector<ClusterType> &dsvToFill,
0089                          edm::RefProd<edmNew::DetSetVector<ClusterType> > &refprod);
0090 
0091     //--- Information about the cloned clusters
0092     std::vector<PixelClusterHitRecord> pixelClusterRecords_;
0093     std::vector<StripClusterHitRecord> stripClusterRecords_;
0094     std::vector<Phase2OTClusterHitRecord> phase2OTClusterRecords_;
0095   };
0096 
0097 }  // namespace helper
0098 
0099 #endif