CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/RecoTracker/TrackProducer/src/ClusterRemovalRefSetter.cc

Go to the documentation of this file.
00001 #include "RecoTracker/TrackProducer/interface/ClusterRemovalRefSetter.h"
00002 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
00003 
00004 ClusterRemovalRefSetter::ClusterRemovalRefSetter(const edm::Event &iEvent, const edm::InputTag tag) {
00005     edm::Handle<reco::ClusterRemovalInfo> hCRI;
00006     iEvent.getByLabel(tag, hCRI);
00007     cri_ = &*hCRI; 
00008 
00009     //std::cout << "Rekeying PIXEL ProdID " << cri_->pixelNewRefProd().id() << " => " << cri_->pixelRefProd().id() << std::endl;
00010     //std::cout << "Rekeying STRIP ProdID " << cri_->stripNewRefProd().id() << " => " << cri_->stripRefProd().id() << std::endl;
00011 }
00012 
00013 void ClusterRemovalRefSetter::reKey(TrackingRecHit *hit) const {
00014     if (!hit->isValid()) return;
00015     DetId detid = hit->geographicalId(); 
00016     uint32_t subdet = detid.subdetId();
00017     if ((subdet == PixelSubdetector::PixelBarrel) || (subdet == PixelSubdetector::PixelEndcap)) {
00018         if (!cri_->hasPixel()) return;
00019         reKeyPixel(reinterpret_cast<SiPixelRecHit *>(hit)->omniCluster());
00020     } else {
00021         if (!cri_->hasStrip()) return;
00022         const std::type_info &type = typeid(*hit);
00023         if (type == typeid(SiStripRecHit2D)) {
00024             reKeyStrip(reinterpret_cast<SiStripRecHit2D *>(hit)->omniCluster());
00025         } else if (type == typeid(SiStripRecHit1D)) {
00026             reKeyStrip(reinterpret_cast<SiStripRecHit1D *>(hit)->omniCluster());
00027         } else if (type == typeid(SiStripMatchedRecHit2D)) {
00028             SiStripMatchedRecHit2D *mhit = reinterpret_cast<SiStripMatchedRecHit2D *>(hit);
00029             reKeyStrip(mhit->monoClusterRef());
00030             reKeyStrip(mhit->stereoClusterRef());
00031         } else if (type == typeid(ProjectedSiStripRecHit2D)) {
00032             ProjectedSiStripRecHit2D *phit = reinterpret_cast<ProjectedSiStripRecHit2D *>(hit);
00033             reKeyStrip(phit->originalHit().omniCluster());
00034         } else throw cms::Exception("Unknown RecHit Type") << "RecHit of type " << type.name() << " not supported. (use c++filt to demangle the name)";
00035     }
00036 }
00037 
00038 
00039 void ClusterRemovalRefSetter::reKeyPixel(OmniClusterRef& newRef) const {
00040   // "newRef" as it refs to the "new"(=cleaned) collection, instead of the old one
00041   using reco::ClusterRemovalInfo;
00042   const ClusterRemovalInfo::Indices &indices = cri_->pixelIndices();
00043   
00044   if (newRef.id()  != cri_->pixelNewRefProd().id()) {
00045     throw cms::Exception("Inconsistent Data") << "ClusterRemovalRefSetter: " << 
00046       "Existing pixel cluster refers to product ID " << newRef.id() << 
00047       " while the ClusterRemovalInfo expects as *new* cluster collection the ID " << cri_->pixelNewRefProd().id() << "\n";
00048   }
00049   size_t newIndex = newRef.key();
00050   assert(newIndex < indices.size());
00051   size_t oldIndex = indices[newIndex];
00052   ClusterPixelRef oldRef(cri_->pixelRefProd(), oldIndex);
00053   newRef = OmniClusterRef(oldRef);
00054 }
00055 
00056 
00057 void ClusterRemovalRefSetter::reKeyStrip(OmniClusterRef& newRef) const {
00058   // "newRef" as it refs to the "new"(=cleaned) collection, instead of the old one
00059   using reco::ClusterRemovalInfo;
00060   const ClusterRemovalInfo::Indices &indices = cri_->stripIndices();
00061   
00062   if (newRef.id() != cri_->stripNewRefProd().id()) {   // this is a cfg error in the tracking configuration, much more likely
00063     throw cms::Exception("Inconsistent Data") << "ClusterRemovalRefSetter: " << 
00064       "Existing strip cluster refers to product ID " << newRef.id() << 
00065       " while the ClusterRemovalInfo expects as *new* cluster collection the ID " << cri_->stripNewRefProd().id() << "\n";
00066   }
00067   
00068   size_t newIndex = newRef.key();
00069   assert(newIndex < indices.size());
00070   size_t oldIndex = indices[newIndex];
00071   ClusterStripRef oldRef(cri_->stripRefProd(), oldIndex);
00072   newRef = OmniClusterRef(oldRef);
00073 }
00074 
00075