CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/CommonTools/RecoAlgos/src/ClusterStorer.cc

Go to the documentation of this file.
00001 #include "CommonTools/RecoAlgos/interface/ClusterStorer.h"
00002 
00003 #include "FWCore/Utilities/interface/Exception.h"
00004 
00005 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
00006 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
00007 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit1D.h"
00008 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2D.h"
00009 #include "DataFormats/TrackerRecHit2D/interface/ProjectedSiStripRecHit2D.h"
00010 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2D.h"
00011 // FastSim hits:
00012 #include "DataFormats/TrackerRecHit2D/interface/SiTrackerGSRecHit2D.h"
00013 #include "DataFormats/TrackerRecHit2D/interface/SiTrackerGSMatchedRecHit2D.h"
00014 
00015 
00016 #include "DataFormats/SiStripDetId/interface/SiStripDetId.h"
00017 
00018 namespace helper {
00019 
00020   // -------------------------------------------------------------
00021   void ClusterStorer::addCluster(TrackingRecHitCollection &hits, size_t index)
00022   {
00023     TrackingRecHit &newHit = hits[index];
00024     const std::type_info &hit_type = typeid(newHit);
00025     if (hit_type == typeid(SiPixelRecHit)) {
00026       //std::cout << "|  It is a Pixel hit !!" << std::endl;
00027       pixelClusterRecords_.push_back(PixelClusterHitRecord(static_cast<SiPixelRecHit&>(newHit),
00028                                                            hits, index));
00029     } else if (hit_type == typeid(SiStripRecHit1D)) {
00030       //std::cout << "|   It is a SiStripRecHit1D hit !!" << std::endl;
00031       stripClusterRecords_.push_back(StripClusterHitRecord(static_cast<SiStripRecHit1D&>(newHit),
00032                                                            hits, index));
00033     } else if (hit_type == typeid(SiStripRecHit2D)) {
00034       //std::cout << "|   It is a SiStripRecHit2D hit !!" << std::endl;
00035       stripClusterRecords_.push_back(StripClusterHitRecord(static_cast<SiStripRecHit2D&>(newHit),
00036                                                            hits, index));
00037     } else if (hit_type == typeid(SiStripMatchedRecHit2D)) {      
00038       //std::cout << "|   It is a SiStripMatchedRecHit2D hit !!" << std::endl;
00039       SiStripMatchedRecHit2D &mhit = static_cast<SiStripMatchedRecHit2D&>(newHit);
00040       stripClusterRecords_.push_back(StripClusterHitRecord(*mhit.monoHit(), hits, index));
00041       stripClusterRecords_.push_back(StripClusterHitRecord(*mhit.stereoHit(), hits, index));
00042     } else if (hit_type == typeid(ProjectedSiStripRecHit2D)) {
00043       //std::cout << "|   It is a ProjectedSiStripRecHit2D hit !!" << std::endl;
00044       ProjectedSiStripRecHit2D &phit = static_cast<ProjectedSiStripRecHit2D&>(newHit);
00045       stripClusterRecords_.push_back(StripClusterHitRecord(phit.originalHit(), hits, index));
00046     } else {
00047       if (hit_type == typeid(SiTrackerGSMatchedRecHit2D)
00048           || hit_type == typeid(SiTrackerGSRecHit2D)) {
00049         //std::cout << "|   It is a " << hit_type.name() << " hit !!" << std::endl;
00050         // FastSim hits: Do nothing instead of caring about FastSim clusters, 
00051         //               not even sure whether these really exist.
00052         //               At least predecessor code in TrackSelector and MuonSelector
00053         //               did not treat them.
00054       } else {
00055         // through for unknown types
00056         throw cms::Exception("UnknownHitType") << "helper::ClusterStorer::addCluster: "
00057                                                << "Unknown hit type " << hit_type.name()
00058                                                << ".\n";
00059       }
00060     } // end 'switch' on hit type
00061     
00062   }
00063   
00064   // -------------------------------------------------------------
00065   void ClusterStorer::clear()
00066   {
00067     pixelClusterRecords_.clear();
00068     stripClusterRecords_.clear();
00069   }
00070   
00071   // -------------------------------------------------------------
00072   void ClusterStorer::
00073   processAllClusters(edmNew::DetSetVector<SiPixelCluster> &pixelDsvToFill,
00074                      edm::RefProd<edmNew::DetSetVector<SiPixelCluster> > refPixelClusters,
00075                      edmNew::DetSetVector<SiStripCluster> &stripDsvToFill,
00076                      edm::RefProd<edmNew::DetSetVector<SiStripCluster> > refStripClusters)
00077   {
00078     if (!pixelClusterRecords_.empty()) {
00079       this->processClusters<SiPixelRecHit, SiPixelCluster>
00080         (pixelClusterRecords_, pixelDsvToFill, refPixelClusters);
00081     }
00082     if (!stripClusterRecords_.empty()) {
00083       // All we need from the HitType 'SiStripRecHit2D' is the
00084       // typedef for 'SiStripRecHit2D::ClusterRef'.
00085       // The fact that rekey<SiStripRecHit2D> is called is irrelevant since
00086       // ClusterHitRecord<typename SiStripRecHit2D::ClusterRef>::rekey<RecHitType>
00087       // is specialised such that 'RecHitType' is not used...
00088       this->processClusters<SiStripRecHit2D, SiStripCluster>
00089         (stripClusterRecords_, stripDsvToFill, refStripClusters);
00090     }
00091   }
00092   
00093   //-------------------------------------------------------------
00094   template<typename HitType, typename ClusterType>
00095   void ClusterStorer::
00096   processClusters(std::vector<ClusterHitRecord<typename HitType::ClusterRef> > &clusterRecords,
00097                   edmNew::DetSetVector<ClusterType>                            &dsvToFill,
00098                   edm::RefProd< edmNew::DetSetVector<ClusterType> >            &refprod)
00099   {
00100     std::sort(clusterRecords.begin(), clusterRecords.end()); // this sorts them by detid 
00101     typedef
00102       typename std::vector<ClusterHitRecord<typename HitType::ClusterRef> >::const_iterator
00103       RIT;
00104     RIT it = clusterRecords.begin(), end = clusterRecords.end();
00105     size_t clusters = 0;
00106     while (it != end) {
00107       RIT it2 = it;
00108       uint32_t detid = it->detid();
00109       
00110       // first isolate all clusters on the same detid
00111       while ( (it2 != end) && (it2->detid() == detid)) {  ++it2; }
00112       // now [it, it2] bracket one detid
00113       
00114       // then prepare to copy the clusters
00115       typename edmNew::DetSetVector<ClusterType>::FastFiller filler(dsvToFill, detid);
00116       typename HitType::ClusterRef lastRef, newRef;
00117       for ( ; it != it2; ++it) { // loop on the detid
00118         // first check if we need to clone the hit
00119         if (it->clusterRef() != lastRef) { 
00120           lastRef = it->clusterRef();
00121           // clone cluster
00122           filler.push_back( *lastRef );  
00123           // make new ref
00124           newRef = typename HitType::ClusterRef( refprod, clusters++ );
00125         } 
00126         it->template rekey<HitType>(newRef);
00127       } // end of the loop on a single detid
00128       
00129     } // end of the loop on all clusters
00130   }
00131   
00132   //-------------------------------------------------------------
00133   // helper classes  
00134   //-------------------------------------------------------------
00135   // generic rekey (in practise for pixel only...)
00136   template<typename ClusterRefType> // template for class
00137   template<typename RecHitType>     // template for member function
00138   void ClusterStorer::ClusterHitRecord<ClusterRefType>::
00139   rekey(const ClusterRefType &newRef) const
00140   {
00141     TrackingRecHit & genericHit = (*hits_)[index_]; 
00142     RecHitType *hit = 0;
00143     if (genericHit.geographicalId().rawId() == detid_) { // a hit on this det, so it's simple
00144       hit = dynamic_cast<RecHitType *>(&genericHit); //static_cast<RecHitType *>(&genericHit);
00145     }
00146     assert (hit != 0);
00147     assert (hit->cluster() == ref_); // otherwise something went wrong
00148     hit->setClusterRef(newRef);
00149   }
00150   
00151   // -------------------------------------------------------------
00152   // Specific rekey for class template ClusterRefType = SiStripRecHit2D::ClusterRef,
00153   // RecHitType is not used.
00154   template<>
00155   template<typename RecHitType> // or template<> to specialise also here?
00156   void ClusterStorer::ClusterHitRecord<SiStripRecHit2D::ClusterRef>::
00157   //  rekey<SiStripRecHit2D>(const SiStripRecHit2D::ClusterRef &newRef) const
00158   rekey(const SiStripRecHit2D::ClusterRef &newRef) const
00159   {
00160     TrackingRecHit &genericHit = (*hits_)[index_];
00161     const std::type_info &hit_type = typeid(genericHit);
00162 
00163     if (typeid(SiStripRecHit1D) == hit_type) {
00164       // case of SiStripRecHit1D 
00165       SiStripRecHit1D *hit = &static_cast<SiStripRecHit1D&>(genericHit);
00166       assert(hit->cluster() == ref_); // otherwise something went wrong
00167       hit->setClusterRef(newRef);
00168     } else {
00169       // various cases resolving to SiStripRecHit2D 
00170       SiStripRecHit2D *hit = 0;
00171       if (typeid(SiStripRecHit2D) == hit_type) {
00172         hit = &static_cast<SiStripRecHit2D&>(genericHit);
00173       } else if (typeid(SiStripMatchedRecHit2D) == hit_type) {
00174         SiStripMatchedRecHit2D &mhit = static_cast<SiStripMatchedRecHit2D&>(genericHit);
00175         hit = (SiStripDetId(detid_).stereo() ? mhit.stereoHit() : mhit.monoHit());
00176       } else if (typeid(ProjectedSiStripRecHit2D) == hit_type) {
00177         ProjectedSiStripRecHit2D &phit = static_cast<ProjectedSiStripRecHit2D&>(genericHit);
00178         hit = &phit.originalHit();
00179       }
00180       assert(hit != 0); // to catch missing RecHit types
00181       assert(hit->cluster() == ref_); // otherwise something went wrong
00182       hit->setClusterRef(newRef);
00183     }
00184 
00185   }
00186   
00187 } // end namespace 'helper'