Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:31:33

0001 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
0002 
0003 ClusterRemovalRefSetter::ClusterRemovalRefSetter(const edm::Event &iEvent, const edm::InputTag &tag) {
0004   edm::Handle<reco::ClusterRemovalInfo> hCRI;
0005   iEvent.getByLabel(tag, hCRI);
0006   cri_ = &*hCRI;
0007 
0008   //std::cout << "Rekeying PIXEL ProdID " << cri_->pixelNewRefProd().id() << " => " << cri_->pixelRefProd().id() << std::endl;
0009   //std::cout << "Rekeying STRIP ProdID " << cri_->stripNewRefProd().id() << " => " << cri_->stripRefProd().id() << std::endl;
0010 }
0011 
0012 ClusterRemovalRefSetter::ClusterRemovalRefSetter(const edm::Event &iEvent,
0013                                                  const edm::EDGetTokenT<reco::ClusterRemovalInfo> &token) {
0014   edm::Handle<reco::ClusterRemovalInfo> hCRI;
0015   iEvent.getByToken(token, hCRI);
0016   cri_ = &*hCRI;
0017 }
0018 
0019 void ClusterRemovalRefSetter::reKey(TrackingRecHit *hit) const {
0020   if (!hit->isValid())
0021     return;
0022   DetId detid = hit->geographicalId();
0023   uint32_t subdet = detid.subdetId();
0024   if ((subdet == PixelSubdetector::PixelBarrel) || (subdet == PixelSubdetector::PixelEndcap)) {
0025     if (!cri_->hasPixel())
0026       return;
0027     reKeyPixel(reinterpret_cast<SiPixelRecHit *>(hit)->omniCluster());
0028   } else {
0029     if (!cri_->hasStrip())
0030       return;
0031     const std::type_info &type = typeid(*hit);
0032     if (type == typeid(SiStripRecHit2D)) {
0033       reKeyStrip(reinterpret_cast<SiStripRecHit2D *>(hit)->omniCluster());
0034     } else if (type == typeid(SiStripRecHit1D)) {
0035       reKeyStrip(reinterpret_cast<SiStripRecHit1D *>(hit)->omniCluster());
0036     } else if (type == typeid(SiStripMatchedRecHit2D)) {
0037       SiStripMatchedRecHit2D *mhit = reinterpret_cast<SiStripMatchedRecHit2D *>(hit);
0038       reKeyStrip(mhit->monoClusterRef());
0039       reKeyStrip(mhit->stereoClusterRef());
0040     } else if (type == typeid(ProjectedSiStripRecHit2D)) {
0041       ProjectedSiStripRecHit2D *phit = reinterpret_cast<ProjectedSiStripRecHit2D *>(hit);
0042       reKeyStrip(phit->originalHit().omniCluster());
0043     } else
0044       throw cms::Exception("Unknown RecHit Type")
0045           << "RecHit of type " << type.name() << " not supported. (use c++filt to demangle the name)";
0046   }
0047 }
0048 
0049 void ClusterRemovalRefSetter::reKeyPixel(OmniClusterRef &newRef) const {
0050   // "newRef" as it refs to the "new"(=cleaned) collection, instead of the old one
0051   using reco::ClusterRemovalInfo;
0052   const ClusterRemovalInfo::Indices &indices = cri_->pixelIndices();
0053 
0054   if (newRef.id() != cri_->pixelNewRefProd().id()) {
0055     throw cms::Exception("Inconsistent Data")
0056         << "ClusterRemovalRefSetter: "
0057         << "Existing pixel cluster refers to product ID " << newRef.id()
0058         << " while the ClusterRemovalInfo expects as *new* cluster collection the ID " << cri_->pixelNewRefProd().id()
0059         << "\n";
0060   }
0061   size_t newIndex = newRef.key();
0062   assert(newIndex < indices.size());
0063   size_t oldIndex = indices[newIndex];
0064   ClusterPixelRef oldRef(cri_->pixelRefProd(), oldIndex);
0065   newRef = OmniClusterRef(oldRef);
0066 }
0067 
0068 void ClusterRemovalRefSetter::reKeyStrip(OmniClusterRef &newRef) const {
0069   // "newRef" as it refs to the "new"(=cleaned) collection, instead of the old one
0070   using reco::ClusterRemovalInfo;
0071   const ClusterRemovalInfo::Indices &indices = cri_->stripIndices();
0072 
0073   if (newRef.id() !=
0074       cri_->stripNewRefProd().id()) {  // this is a cfg error in the tracking configuration, much more likely
0075     throw cms::Exception("Inconsistent Data")
0076         << "ClusterRemovalRefSetter: "
0077         << "Existing strip cluster refers to product ID " << newRef.id()
0078         << " while the ClusterRemovalInfo expects as *new* cluster collection the ID " << cri_->stripNewRefProd().id()
0079         << "\n";
0080   }
0081 
0082   size_t newIndex = newRef.key();
0083   assert(newIndex < indices.size());
0084   size_t oldIndex = indices[newIndex];
0085   ClusterStripRef oldRef(cri_->stripRefProd(), oldIndex);
0086   newRef = OmniClusterRef(oldRef);
0087 }