Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:01:10

0001 #include "CommonTools/RecoAlgos/interface/ClusterStorer.h"
0002 
0003 #include "FWCore/Utilities/interface/Exception.h"
0004 
0005 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
0006 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
0007 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit1D.h"
0008 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2D.h"
0009 #include "DataFormats/TrackerRecHit2D/interface/ProjectedSiStripRecHit2D.h"
0010 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2D.h"
0011 #include "DataFormats/TrackerRecHit2D/interface/Phase2TrackerRecHit1D.h"
0012 #include "DataFormats/TrackerRecHit2D/interface/VectorHit.h"
0013 // FastSim hits:
0014 #include "DataFormats/TrackerRecHit2D/interface/FastTrackerRecHit.h"
0015 #include "DataFormats/TrackerRecHit2D/interface/FastProjectedTrackerRecHit.h"
0016 #include "DataFormats/TrackerRecHit2D/interface/FastMatchedTrackerRecHit.h"
0017 
0018 #include "DataFormats/SiStripDetId/interface/SiStripDetId.h"
0019 
0020 namespace helper {
0021 
0022   // -------------------------------------------------------------
0023   //FIXME (push cluster pointers...)
0024   void ClusterStorer::addCluster(TrackingRecHitCollection &hits, size_t index) {
0025     TrackingRecHit &newHit = hits[index];
0026     const std::type_info &hit_type = typeid(newHit);
0027     if (hit_type == typeid(SiPixelRecHit)) {
0028       //std::cout << "|  It is a Pixel hit !!" << std::endl;
0029       pixelClusterRecords_.push_back(PixelClusterHitRecord(static_cast<SiPixelRecHit &>(newHit), hits, index));
0030     } else if (hit_type == typeid(SiStripRecHit1D)) {
0031       //std::cout << "|   It is a SiStripRecHit1D hit !!" << std::endl;
0032       stripClusterRecords_.push_back(StripClusterHitRecord(static_cast<SiStripRecHit1D &>(newHit), hits, index));
0033     } else if (hit_type == typeid(SiStripRecHit2D)) {
0034       //std::cout << "|   It is a SiStripRecHit2D hit !!" << std::endl;
0035       stripClusterRecords_.push_back(StripClusterHitRecord(static_cast<SiStripRecHit2D &>(newHit), hits, index));
0036     } else if (hit_type == typeid(SiStripMatchedRecHit2D)) {
0037       //std::cout << "|   It is a SiStripMatchedRecHit2D hit !!" << std::endl;
0038       SiStripMatchedRecHit2D &mhit = static_cast<SiStripMatchedRecHit2D &>(newHit);
0039       stripClusterRecords_.push_back(StripClusterHitRecord(mhit.monoHit(), hits, index));
0040       stripClusterRecords_.push_back(StripClusterHitRecord(mhit.stereoHit(), hits, index));
0041     } else if (hit_type == typeid(ProjectedSiStripRecHit2D)) {
0042       //std::cout << "|   It is a ProjectedSiStripRecHit2D hit !!" << std::endl;
0043       ProjectedSiStripRecHit2D &phit = static_cast<ProjectedSiStripRecHit2D &>(newHit);
0044       stripClusterRecords_.push_back(StripClusterHitRecord(phit.originalHit(), hits, index));
0045     } else if (hit_type == typeid(Phase2TrackerRecHit1D)) {
0046       //FIXME:: this is just temporary solution for phase2,
0047       //it is not really running in the phase2 tracking wf - yet...
0048       phase2OTClusterRecords_.push_back(
0049           Phase2OTClusterHitRecord(static_cast<Phase2TrackerRecHit1D &>(newHit), hits, index));
0050     } else if (hit_type == typeid(VectorHit)) {
0051       //FIXME:: this is just temporary solution for phase2,
0052       //the VectorHit has 2 clusters but just a hit!
0053       phase2OTClusterRecords_.push_back(Phase2OTClusterHitRecord(static_cast<VectorHit &>(newHit), hits, index));
0054     } else {
0055       if (hit_type == typeid(FastTrackerRecHit) || hit_type == typeid(FastProjectedTrackerRecHit) ||
0056           hit_type == typeid(FastMatchedTrackerRecHit)) {
0057         //std::cout << "|   It is a " << hit_type.name() << " hit !!" << std::endl;
0058         // FastSim hits: Do nothing instead of caring about FastSim clusters,
0059         //               not even sure whether these really exist.
0060         //               At least predecessor code in TrackSelector and MuonSelector
0061         //               did not treat them.
0062       } else {
0063         // through for unknown types
0064         throw cms::Exception("UnknownHitType") << "helper::ClusterStorer::addCluster: "
0065                                                << "Unknown hit type " << hit_type.name() << ".\n";
0066       }
0067     }  // end 'switch' on hit type
0068   }
0069 
0070   // -------------------------------------------------------------
0071   void ClusterStorer::clear() {
0072     pixelClusterRecords_.clear();
0073     stripClusterRecords_.clear();
0074     phase2OTClusterRecords_.clear();
0075   }
0076 
0077   // -------------------------------------------------------------
0078   void ClusterStorer::processAllClusters(
0079       edmNew::DetSetVector<SiPixelCluster> &pixelDsvToFill,
0080       edm::RefProd<edmNew::DetSetVector<SiPixelCluster> > refPixelClusters,
0081       edmNew::DetSetVector<SiStripCluster> &stripDsvToFill,
0082       edm::RefProd<edmNew::DetSetVector<SiStripCluster> > refStripClusters,
0083       edmNew::DetSetVector<Phase2TrackerCluster1D> &phase2OTDsvToFill,
0084       edm::RefProd<edmNew::DetSetVector<Phase2TrackerCluster1D> > refPhase2OTClusters) {
0085     if (!pixelClusterRecords_.empty()) {
0086       this->processClusters<SiPixelRecHit, SiPixelCluster>(pixelClusterRecords_, pixelDsvToFill, refPixelClusters);
0087     }
0088     if (!stripClusterRecords_.empty()) {
0089       // All we need from the HitType 'SiStripRecHit2D' is the
0090       // typedef for 'SiStripRecHit2D::ClusterRef'.
0091       // The fact that rekey<SiStripRecHit2D> is called is irrelevant since
0092       // ClusterHitRecord<typename SiStripRecHit2D::ClusterRef>::rekey<RecHitType>
0093       // is specialised such that 'RecHitType' is not used...
0094       this->processClusters<SiStripRecHit2D, SiStripCluster>(stripClusterRecords_, stripDsvToFill, refStripClusters);
0095     }
0096     if (!phase2OTClusterRecords_.empty()) {
0097       this->processClusters<Phase2TrackerRecHit1D, Phase2TrackerCluster1D>(
0098           phase2OTClusterRecords_, phase2OTDsvToFill, refPhase2OTClusters);
0099     }
0100   }
0101 
0102   //-------------------------------------------------------------
0103   template <typename HitType, typename ClusterType>
0104   void ClusterStorer::processClusters(std::vector<ClusterHitRecord<typename HitType::ClusterRef> > &clusterRecords,
0105                                       edmNew::DetSetVector<ClusterType> &dsvToFill,
0106                                       edm::RefProd<edmNew::DetSetVector<ClusterType> > &refprod) {
0107     std::sort(clusterRecords.begin(), clusterRecords.end());  // this sorts them by detid
0108     typedef typename std::vector<ClusterHitRecord<typename HitType::ClusterRef> >::iterator RIT;
0109     RIT it = clusterRecords.begin(), end = clusterRecords.end();
0110     size_t clusters = 0;
0111     while (it != end) {
0112       RIT it2 = it;
0113       uint32_t detid = it->detid();
0114 
0115       // first isolate all clusters on the same detid
0116       while ((it2 != end) && (it2->detid() == detid)) {
0117         ++it2;
0118       }
0119       // now [it, it2] bracket one detid
0120 
0121       // then prepare to copy the clusters
0122       typename edmNew::DetSetVector<ClusterType>::FastFiller filler(dsvToFill, detid);
0123       typename HitType::ClusterRef lastRef, newRef;
0124       for (; it != it2; ++it) {  // loop on the detid
0125         // first check if we need to clone the hit
0126         if (it->clusterRef() != lastRef) {
0127           lastRef = it->clusterRef();
0128           // clone cluster
0129           filler.push_back(*lastRef);
0130           // make new ref
0131           newRef = typename HitType::ClusterRef(refprod, clusters++);
0132         }
0133         it->template rekey<HitType>(newRef);
0134       }  // end of the loop on a single detid
0135 
0136     }  // end of the loop on all clusters
0137   }
0138 
0139   //-------------------------------------------------------------
0140   // helper classes
0141   //-------------------------------------------------------------
0142   // FIXME (migrate to new RTTI and interface)
0143   // generic rekey (in practise for pixel only...)
0144   template <typename ClusterRefType>  // template for class
0145   template <typename RecHitType>      // template for member function
0146   void ClusterStorer::ClusterHitRecord<ClusterRefType>::rekey(const ClusterRefType &newRef) {
0147     TrackingRecHit &genericHit = (*hits_)[index_];
0148     RecHitType *hit = nullptr;
0149     if (genericHit.geographicalId().rawId() == detid_) {  // a hit on this det, so it's simple
0150       hit = dynamic_cast<RecHitType *>(&genericHit);      //static_cast<RecHitType *>(&genericHit);
0151     }
0152     assert(hit != nullptr);
0153     assert(hit->cluster() == ref_);  // otherwise something went wrong
0154     hit->setClusterRef(newRef);
0155   }
0156 
0157   // -------------------------------------------------------------
0158   // Specific rekey for class template ClusterRefType = SiStripRecHit2D::ClusterRef,
0159   // RecHitType is not used.
0160   template <>
0161   template <typename RecHitType>  // or template<> to specialise also here?
0162   void ClusterStorer::ClusterHitRecord<SiStripRecHit2D::ClusterRef>::rekey(const SiStripRecHit2D::ClusterRef &newRef) {
0163     TrackingRecHit &genericHit = (*hits_)[index_];
0164     const std::type_info &hit_type = typeid(genericHit);
0165 
0166     OmniClusterRef *cluRef = nullptr;
0167     if (typeid(SiStripRecHit1D) == hit_type) {
0168       cluRef = &static_cast<SiStripRecHit1D &>(genericHit).omniCluster();
0169     } else if (typeid(SiStripRecHit2D) == hit_type) {
0170       cluRef = &static_cast<SiStripRecHit2D &>(genericHit).omniCluster();
0171     } else if (typeid(SiStripMatchedRecHit2D) == hit_type) {
0172       SiStripMatchedRecHit2D &mhit = static_cast<SiStripMatchedRecHit2D &>(genericHit);
0173       cluRef = (SiStripDetId(detid_).stereo() ? &mhit.stereoClusterRef() : &mhit.monoClusterRef());
0174     } else if (typeid(ProjectedSiStripRecHit2D) == hit_type) {
0175       cluRef = &static_cast<ProjectedSiStripRecHit2D &>(genericHit).omniCluster();
0176     }
0177 
0178     assert(cluRef != nullptr);            // to catch missing RecHit types
0179     assert(cluRef->key() == ref_.key());  // otherwise something went wrong
0180     (*cluRef) = OmniClusterRef(newRef);
0181   }
0182 
0183   // -------------------------------------------------------------
0184   // Specific rekey for class template ClusterRefType = Phase2TrackerRecHit1D::ClusterRef
0185   template <>
0186   template <typename RecHitType>
0187   void ClusterStorer::ClusterHitRecord<Phase2TrackerRecHit1D::ClusterRef>::rekey(
0188       const Phase2TrackerRecHit1D::ClusterRef &newRef) {
0189     TrackingRecHit &genericHit = (*hits_)[index_];
0190     const std::type_info &hit_type = typeid(genericHit);
0191 
0192     OmniClusterRef *cluRef = nullptr;
0193     if (typeid(Phase2TrackerRecHit1D) == hit_type) {
0194       cluRef = &static_cast<Phase2TrackerRecHit1D &>(genericHit).omniCluster();
0195     } else if (typeid(VectorHit) == hit_type) {
0196       VectorHit &vHit = static_cast<VectorHit &>(genericHit);
0197       // FIXME: this essentially uses a hack
0198       // https://github.com/cms-sw/cmssw/blob/master/DataFormats/TrackerCommon/interface/TrackerTopology.h#L291
0199       cluRef = (SiStripDetId(detid_).stereo() ? &vHit.upperClusterRef() : &vHit.lowerClusterRef());
0200     }
0201 
0202     assert(cluRef != nullptr);            // to catch missing RecHit types
0203     assert(cluRef->key() == ref_.key());  // otherwise something went wrong
0204     (*cluRef) = OmniClusterRef(newRef);
0205   }
0206 
0207 }  // namespace helper