CMS 3D CMS Logo

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     iEvent.get(cri_->pixelProdID(), handlePixel_);
00010     iEvent.get(cri_->stripProdID(), handleStrip_);
00011 
00012     //std::cout << "Rekeying PIXEL ProdID " << cri_->pixelNewProdID() << " => " << cri_->pixelProdID() << std::endl;
00013     //std::cout << "Rekeying STRIP ProdID " << cri_->stripNewProdID() << " => " << cri_->stripProdID() << std::endl;
00014 }
00015 
00016 void ClusterRemovalRefSetter::reKey(TrackingRecHit *hit) const {
00017     if (!hit->isValid()) return;
00018     DetId detid = hit->geographicalId(); 
00019     uint32_t subdet = detid.subdetId();
00020     if ((subdet == PixelSubdetector::PixelBarrel) || (subdet == PixelSubdetector::PixelEndcap)) {
00021         reKey(reinterpret_cast<SiPixelRecHit *>(hit), detid.rawId());
00022     } else {
00023         const std::type_info &type = typeid(*hit);
00024         if (type == typeid(SiStripRecHit2D)) {
00025             reKey(reinterpret_cast<SiStripRecHit2D *>(hit), detid.rawId());
00026         } else if (type == typeid(SiStripMatchedRecHit2D)) {
00027             SiStripMatchedRecHit2D *mhit = reinterpret_cast<SiStripMatchedRecHit2D *>(hit);
00028             // const_cast is needed: monoHit() e stereoHit() are const only - at least for now
00029             reKey(mhit->monoHit(), mhit->monoHit()->geographicalId().rawId());
00030             reKey(mhit->stereoHit(), mhit->stereoHit()->geographicalId().rawId());
00031         } else if (type == typeid(ProjectedSiStripRecHit2D)) {
00032             ProjectedSiStripRecHit2D *phit = reinterpret_cast<ProjectedSiStripRecHit2D *>(hit);
00033             reKey(&phit->originalHit(), phit->originalHit().geographicalId().rawId());
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 void ClusterRemovalRefSetter::reKey(SiStripRecHit2D *hit, uint32_t detid) const {
00039     using reco::ClusterRemovalInfo;
00040     const ClusterRemovalInfo::Indices &indices = cri_->stripIndices();
00041     SiStripRecHit2D::ClusterRef newRef = hit->cluster();
00042     // "newRef" as it refs to the "new"(=cleaned) collection, instead of the old one
00043 
00044     if (newRef.id() != cri_->stripNewProdID()) {   // this is a cfg error in the tracking configuration, much more likely
00045         throw cms::Exception("Inconsistent Data") << "ClusterRemovalRefSetter: " << 
00046             "Existing strip cluster refers to product ID " << newRef.id() << 
00047             " while the ClusterRemovalInfo expects as *new* cluster collection the ID " << cri_->stripNewProdID() << "\n";
00048     }
00049 
00050     size_t newIndex = newRef.key();
00051     assert(newIndex < indices.size());
00052     size_t oldIndex = indices[newIndex];
00053     SiStripRecHit2D::ClusterRef oldRef(handleStrip_, oldIndex, false);
00054     hit->setClusterRef(oldRef);
00055 }
00056 
00057 void ClusterRemovalRefSetter::reKey(SiPixelRecHit *hit, uint32_t detid) const {
00058     using reco::ClusterRemovalInfo;
00059     const ClusterRemovalInfo::Indices &indices = cri_->pixelIndices();
00060     SiPixelRecHit::ClusterRef newRef = hit->cluster();
00061     // "newRef" as it refs to the "new"(=cleaned) collection, instead of the old one
00062 
00063     if (newRef.id()  != cri_->pixelNewProdID()) {
00064         throw cms::Exception("Inconsistent Data") << "ClusterRemovalRefSetter: " << 
00065             "Existing pixel cluster refers to product ID " << newRef.id() << 
00066             " while the ClusterRemovalInfo expects as *new* cluster collection the ID " << cri_->pixelNewProdID() << "\n";
00067     }
00068     size_t newIndex = newRef.key();
00069     assert(newIndex < indices.size());
00070     size_t oldIndex = indices[newIndex];
00071     SiPixelRecHit::ClusterRef oldRef(handlePixel_, oldIndex, false);
00072     hit->setClusterRef(oldRef);
00073 }
00074 
00075 

Generated on Tue Jun 9 17:46:00 2009 for CMSSW by  doxygen 1.5.4