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
00010
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 reKey(reinterpret_cast<SiPixelRecHit *>(hit), detid.rawId());
00020 } else {
00021 if (!cri_->hasStrip()) return;
00022 const std::type_info &type = typeid(*hit);
00023 if (type == typeid(SiStripRecHit2D)) {
00024 reKey(reinterpret_cast<SiStripRecHit2D *>(hit), detid.rawId());
00025 } else if (type == typeid(SiStripRecHit1D)) {
00026 reKey(reinterpret_cast<SiStripRecHit1D *>(hit), detid.rawId());
00027 } else if (type == typeid(SiStripMatchedRecHit2D)) {
00028 SiStripMatchedRecHit2D *mhit = reinterpret_cast<SiStripMatchedRecHit2D *>(hit);
00029
00030 reKey(mhit->monoHit(), mhit->monoHit()->geographicalId().rawId());
00031 reKey(mhit->stereoHit(), mhit->stereoHit()->geographicalId().rawId());
00032 } else if (type == typeid(ProjectedSiStripRecHit2D)) {
00033 ProjectedSiStripRecHit2D *phit = reinterpret_cast<ProjectedSiStripRecHit2D *>(hit);
00034 reKey(&phit->originalHit(), phit->originalHit().geographicalId().rawId());
00035 } else throw cms::Exception("Unknown RecHit Type") << "RecHit of type " << type.name() << " not supported. (use c++filt to demangle the name)";
00036 }
00037 }
00038
00039 void ClusterRemovalRefSetter::reKey(SiStripRecHit2D *hit, uint32_t detid) const {
00040 using reco::ClusterRemovalInfo;
00041 const ClusterRemovalInfo::Indices &indices = cri_->stripIndices();
00042 SiStripRecHit2D::ClusterRef newRef = hit->cluster();
00043
00044
00045 if (newRef.id() != cri_->stripNewRefProd().id()) {
00046 throw cms::Exception("Inconsistent Data") << "ClusterRemovalRefSetter: " <<
00047 "Existing strip cluster refers to product ID " << newRef.id() <<
00048 " while the ClusterRemovalInfo expects as *new* cluster collection the ID " << cri_->stripNewRefProd().id() << "\n";
00049 }
00050
00051 size_t newIndex = newRef.key();
00052 assert(newIndex < indices.size());
00053 size_t oldIndex = indices[newIndex];
00054 SiStripRecHit2D::ClusterRef oldRef(cri_->stripRefProd(), oldIndex);
00055 hit->setClusterRef(oldRef);
00056 }
00057
00058 void ClusterRemovalRefSetter::reKey(SiStripRecHit1D *hit, uint32_t detid) const {
00059 using reco::ClusterRemovalInfo;
00060 const ClusterRemovalInfo::Indices &indices = cri_->stripIndices();
00061 SiStripRecHit1D::ClusterRef newRef = hit->cluster();
00062
00063
00064 if (newRef.id() != cri_->stripNewRefProd().id()) {
00065 throw cms::Exception("Inconsistent Data") << "ClusterRemovalRefSetter: " <<
00066 "Existing strip cluster refers to product ID " << newRef.id() <<
00067 " while the ClusterRemovalInfo expects as *new* cluster collection the ID " << cri_->stripNewRefProd().id() << "\n";
00068 }
00069
00070 size_t newIndex = newRef.key();
00071 assert(newIndex < indices.size());
00072 size_t oldIndex = indices[newIndex];
00073 SiStripRecHit1D::ClusterRef oldRef(cri_->stripRefProd(), oldIndex);
00074 hit->setClusterRef(oldRef);
00075 }
00076
00077 void ClusterRemovalRefSetter::reKey(SiPixelRecHit *hit, uint32_t detid) const {
00078 using reco::ClusterRemovalInfo;
00079 const ClusterRemovalInfo::Indices &indices = cri_->pixelIndices();
00080 SiPixelRecHit::ClusterRef newRef = hit->cluster();
00081
00082
00083 if (newRef.id() != cri_->pixelNewRefProd().id()) {
00084 throw cms::Exception("Inconsistent Data") << "ClusterRemovalRefSetter: " <<
00085 "Existing pixel cluster refers to product ID " << newRef.id() <<
00086 " while the ClusterRemovalInfo expects as *new* cluster collection the ID " << cri_->pixelNewRefProd().id() << "\n";
00087 }
00088 size_t newIndex = newRef.key();
00089 assert(newIndex < indices.size());
00090 size_t oldIndex = indices[newIndex];
00091 SiPixelRecHit::ClusterRef oldRef(cri_->pixelRefProd(), oldIndex);
00092 hit->setClusterRef(oldRef);
00093 }
00094
00095