CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/RecoTracker/TkSeedingLayers/src/HitExtractorSTRP.cc

Go to the documentation of this file.
00001 #include "HitExtractorSTRP.h"
00002 #include "Geometry/TrackerGeometryBuilder/interface/TrackerLayerIdAccessor.h"
00003 #include "TrackingTools/DetLayers/interface/DetLayer.h"
00004 
00005 #include "DataFormats/Common/interface/Handle.h"
00006 #include "FWCore/Framework/interface/Event.h"
00007 
00008 #include "FWCore/Framework/interface/EventSetup.h"
00009 #include "FWCore/Framework/interface/ESHandle.h"
00010 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2DCollection.h"
00011 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2DCollection.h"
00012 
00013 #include "DataFormats/SiStripDetId/interface/TIBDetId.h"
00014 #include "DataFormats/SiStripDetId/interface/TIDDetId.h"
00015 #include "DataFormats/SiStripDetId/interface/TOBDetId.h"
00016 #include "DataFormats/SiStripDetId/interface/TECDetId.h"
00017 
00018 #include "TrackingTools/TransientTrackingRecHit/interface/TrackingRecHitProjector.h"
00019 #include "RecoTracker/TransientTrackingRecHit/interface/ProjectedRecHit2D.h"
00020 
00021 using namespace ctfseeding;
00022 using namespace std;
00023 using namespace edm;
00024 
00025 HitExtractorSTRP::HitExtractorSTRP( const DetLayer* detLayer, 
00026     SeedingLayer::Side & side, int idLayer)
00027   : theLayer(detLayer), theSide(side), theIdLayer(idLayer),
00028     hasMatchedHits(false), hasRPhiHits(false), hasStereoHits(false),
00029     hasRingSelector(false), theMinRing(1), theMaxRing(0), hasSimpleRphiHitsCleaner(true)
00030 { }
00031 
00032 void HitExtractorSTRP::useRingSelector(int minRing, int maxRing) 
00033 {
00034   hasRingSelector=true;
00035   theMinRing=minRing;
00036   theMaxRing=maxRing; 
00037 }
00038 
00039 bool HitExtractorSTRP::ringRange(int ring) const
00040 {
00041   if (!hasRingSelector) return true;
00042   else if ( ring >= theMinRing && ring <= theMaxRing) return true;
00043   else return false;
00044 }
00045 
00046 bool HitExtractorSTRP::skipThis(const SiStripRecHit2D * hit,
00047                                 edm::Handle<edmNew::DetSetVector<SiStripClusterRef> > & stripClusterRefs) const {
00048   static DetId lastId=hit->geographicalId();
00049   static edmNew::DetSetVector<SiStripClusterRef>::const_iterator f=stripClusterRefs->find(lastId.rawId());
00050   if (hit->geographicalId()!=lastId){
00051     lastId=hit->geographicalId();
00052     f=stripClusterRefs->find(lastId.rawId());
00053   }  
00054   if (f==stripClusterRefs->end()) return false;
00055   if (!hit->isValid())  return false;
00056 
00057   bool skipping=(find(f->begin(),f->end(),hit->cluster())!=f->end());
00058   //if (skipping) LogDebug("HitExtractorSTRP")<<"skipping a hit on :"<<hit->geographicalId().rawId()<<" key: "<<hit->cluster().key();
00059   return skipping;
00060 }
00061 
00062 
00063 void HitExtractorSTRP::project(TransientTrackingRecHit::ConstRecHitPointer & ptr,
00064                                const SiStripRecHit2D * hit,
00065                                TransientTrackingRecHit::ConstRecHitPointer & replaceMe) const{
00066   
00067   if (failProjection) {replaceMe=0; return;}
00068   TrackingRecHitProjector<ProjectedRecHit2D> proj;
00069   TransientTrackingRecHit::RecHitPointer sHit=theSLayer->hitBuilder()->build(hit);
00070   replaceMe=proj.project( *sHit, *ptr->det());
00071   if (!replaceMe) LogDebug("HitExtractorSTRP")<<"projection failed.";
00072 }
00073 
00074 bool HitExtractorSTRP::skipThis(TransientTrackingRecHit::ConstRecHitPointer & ptr,
00075                                 edm::Handle<edmNew::DetSetVector<SiStripClusterRef> > & stripClusterRefs,
00076                                 TransientTrackingRecHit::ConstRecHitPointer & replaceMe) const {
00077   const SiStripMatchedRecHit2D * hit = (SiStripMatchedRecHit2D *) ptr->hit();
00078 
00079   bool rejectSt=false,rejectMono=false;
00080   if (skipThis(hit->stereoHit(),stripClusterRefs))  rejectSt=true;
00081   if (skipThis(hit->monoHit(),stripClusterRefs))    rejectMono=true;
00082 
00083   if (rejectSt&&rejectMono){
00084     //only skip if both hits are done
00085     return true;
00086   }
00087   else{
00088     if (rejectSt) project(ptr,hit->stereoHit(),replaceMe);
00089     else if (rejectMono) project(ptr,hit->monoHit(),replaceMe);
00090     if (!replaceMe) return true; //means that the projection failed, and needs to be skipped
00091     if (rejectSt)
00092       LogDebug("HitExtractorSTRP")<<"a matched hit is partially masked, and the mono hit got projected onto: "<<replaceMe->hit()->geographicalId().rawId()<<" key: "<<hit->monoHit()->cluster().key();
00093     else if (rejectMono)
00094       LogDebug("HitExtractorSTRP")<<"a matched hit is partially masked, and the stereo hit got projected onto: "<<replaceMe->hit()->geographicalId().rawId()<<" key: "<<hit->stereoHit()->cluster().key();
00095     return false; //means the projection succeeded or nothing to be masked, no need to skip and replaceMe is going to be used anyways.
00096   }
00097   return false;
00098 }
00099 
00100 
00101 void HitExtractorSTRP::cleanedOfClusters( const edm::Event& ev, HitExtractor::Hits & hits,
00102                                           bool matched)const{
00103   LogDebug("HitExtractorPIX")<<"getting: "<<hits.size()<<" in input.";
00104   edm::Handle<edmNew::DetSetVector<SiStripClusterRef> > stripClusterRefs;
00105   ev.getByLabel(theSkipClusters,stripClusterRefs);
00106   HitExtractor::Hits newHits;
00107   uint skipped=0;
00108   uint projected=0;
00109   newHits.reserve(hits.size());
00110   TransientTrackingRecHit::ConstRecHitPointer replaceMe;
00111   for (unsigned int iH=0;iH!=hits.size();++iH){
00112     replaceMe=hits[iH];
00113     if (matched && skipThis(hits[iH],stripClusterRefs,replaceMe)){
00114       LogDebug("HitExtractorSTRP")<<"skipping a matched hit on :"<<hits[iH]->hit()->geographicalId().rawId();
00115       skipped++;
00116       continue;
00117     }
00118     if (!matched && skipThis((SiStripRecHit2D*) hits[iH]->hit(),stripClusterRefs)){
00119       LogDebug("HitExtractorSTRP")<<"skipping a hit on :"<<hits[iH]->hit()->geographicalId().rawId()<<" key: ";
00120       skipped++;
00121       continue;
00122     }
00123     if (replaceMe!=hits[iH]) projected++;
00124     newHits.push_back(replaceMe);
00125   }
00126   LogDebug("HitExtractorPIX")<<"skipped :"<<skipped<<" strip rechits because of clusters and projected: "<<projected;
00127   hits.swap(newHits);
00128 }
00129 
00130 HitExtractor::Hits HitExtractorSTRP::hits(const SeedingLayer & sl, const edm::Event& ev, const edm::EventSetup& es) const
00131 {
00132   HitExtractor::Hits result;
00133   TrackerLayerIdAccessor accessor;
00134   theSLayer=&sl;
00135   //
00136   // TIB
00137   //
00138   if (theLayer->subDetector() == GeomDetEnumerators::TIB) {
00139     if (hasMatchedHits) {
00140       edm::Handle<SiStripMatchedRecHit2DCollection> matchedHits;
00141       ev.getByLabel( theMatchedHits, matchedHits);
00142       range2SeedingHits( *matchedHits, result, accessor.stripTIBLayer(theIdLayer), sl, es); 
00143       if (skipClusters) cleanedOfClusters(ev,result,true);
00144     }
00145     if (hasRPhiHits) {
00146       edm::Handle<SiStripRecHit2DCollection> rphiHits;
00147       ev.getByLabel( theRPhiHits, rphiHits);
00148       if (hasMatchedHits){ 
00149         if (!hasSimpleRphiHitsCleaner){ // this is a brutal "cleaning". Add something smarter in the future
00150           range2SeedingHits( *rphiHits, result, accessor.stripTIBLayer(theIdLayer), sl, es); 
00151           if (skipClusters) cleanedOfClusters(ev,result,false);
00152         }
00153       } else {
00154         range2SeedingHits( *rphiHits, result, accessor.stripTIBLayer(theIdLayer), sl, es); 
00155         if (skipClusters) cleanedOfClusters(ev,result,false);
00156       }
00157     }
00158     if (hasStereoHits) {
00159       edm::Handle<SiStripRecHit2DCollection> stereoHits;
00160       ev.getByLabel( theStereoHits, stereoHits);
00161       range2SeedingHits( *stereoHits, result, accessor.stripTIBLayer(theIdLayer), sl, es); 
00162       if (skipClusters) cleanedOfClusters(ev,result,false);
00163     }
00164   }
00165   
00166   //
00167   // TID
00168   //
00169   else if (theLayer->subDetector() == GeomDetEnumerators::TID) {
00170       if (hasMatchedHits) {
00171           edm::Handle<SiStripMatchedRecHit2DCollection> matchedHits;
00172           ev.getByLabel( theMatchedHits, matchedHits);
00173           std::pair<DetId,DetIdTIDSameDiskComparator> getter = accessor.stripTIDDisk(theSide,theIdLayer);
00174           SiStripMatchedRecHit2DCollection::Range range = matchedHits->equal_range(getter.first, getter.second);
00175           for (SiStripMatchedRecHit2DCollection::const_iterator it = range.first; it != range.second; ++it) {
00176               int ring = TIDDetId( it->detId() ).ring();  if (!ringRange(ring)) continue;
00177               for (SiStripMatchedRecHit2DCollection::DetSet::const_iterator hit = it->begin(), end = it->end(); hit != end; ++hit) {
00178                 result.push_back( sl.hitBuilder()->build(hit) ); 
00179               }
00180           }
00181           if (skipClusters) cleanedOfClusters(ev,result,true);
00182       }
00183       if (hasRPhiHits) {
00184           edm::Handle<SiStripRecHit2DCollection> rphiHits;
00185           ev.getByLabel( theRPhiHits, rphiHits);
00186           std::pair<DetId,DetIdTIDSameDiskComparator> getter = accessor.stripTIDDisk(theSide,theIdLayer);
00187           SiStripRecHit2DCollection::Range range = rphiHits->equal_range(getter.first, getter.second);
00188           for (SiStripRecHit2DCollection::const_iterator it = range.first; it != range.second; ++it) {
00189               int ring = TIDDetId( it->detId() ).ring();  if (!ringRange(ring)) continue;
00190               if ((SiStripDetId(it->detId()).partnerDetId() != 0) && hasSimpleRphiHitsCleaner) continue;  // this is a brutal "cleaning". Add something smarter in the future
00191               for (SiStripRecHit2DCollection::DetSet::const_iterator hit = it->begin(), end = it->end(); hit != end; ++hit) {
00192                   result.push_back( sl.hitBuilder()->build(hit) );
00193               }
00194           }
00195           if (skipClusters) cleanedOfClusters(ev,result,false);
00196       }
00197       if (hasStereoHits) {
00198           edm::Handle<SiStripRecHit2DCollection> stereoHits;
00199           ev.getByLabel( theStereoHits, stereoHits);
00200           std::pair<DetId,DetIdTIDSameDiskComparator> getter = accessor.stripTIDDisk(theSide,theIdLayer);
00201           SiStripRecHit2DCollection::Range range = stereoHits->equal_range(getter.first, getter.second);
00202           for (SiStripRecHit2DCollection::const_iterator it = range.first; it != range.second; ++it) {
00203               int ring = TIDDetId( it->detId() ).ring();  if (!ringRange(ring)) continue;
00204               for (SiStripRecHit2DCollection::DetSet::const_iterator hit = it->begin(), end = it->end(); hit != end; ++hit) {
00205                   result.push_back( sl.hitBuilder()->build(hit) );
00206               }
00207           }
00208           if (skipClusters) cleanedOfClusters(ev,result,false);
00209       }
00210   }
00211   //
00212   // TOB
00213   //
00214   else if (theLayer->subDetector() == GeomDetEnumerators::TOB) {
00215     if (hasMatchedHits) {
00216       edm::Handle<SiStripMatchedRecHit2DCollection> matchedHits;
00217       ev.getByLabel( theMatchedHits, matchedHits);
00218       range2SeedingHits( *matchedHits, result, accessor.stripTOBLayer(theIdLayer), sl, es); 
00219       if (skipClusters) cleanedOfClusters(ev,result,true);
00220     }
00221     if (hasRPhiHits) {
00222       edm::Handle<SiStripRecHit2DCollection> rphiHits;
00223       ev.getByLabel( theRPhiHits, rphiHits);
00224       if (hasMatchedHits){ 
00225         if (!hasSimpleRphiHitsCleaner){ // this is a brutal "cleaning". Add something smarter in the future
00226           range2SeedingHits( *rphiHits, result, accessor.stripTOBLayer(theIdLayer), sl, es); 
00227           if (skipClusters) cleanedOfClusters(ev,result,false);
00228         }
00229       } else {
00230         range2SeedingHits( *rphiHits, result, accessor.stripTOBLayer(theIdLayer), sl, es); 
00231         if (skipClusters) cleanedOfClusters(ev,result,false);
00232       }
00233     }
00234     if (hasStereoHits) {
00235       edm::Handle<SiStripRecHit2DCollection> stereoHits;
00236       ev.getByLabel( theStereoHits, stereoHits);
00237       range2SeedingHits( *stereoHits, result, accessor.stripTOBLayer(theIdLayer), sl, es); 
00238       if (skipClusters) cleanedOfClusters(ev,result,false);
00239     }
00240   }
00241 
00242   //
00243   // TEC
00244   //
00245   else if (theLayer->subDetector() == GeomDetEnumerators::TEC) {
00246       if (hasMatchedHits) {
00247           edm::Handle<SiStripMatchedRecHit2DCollection> matchedHits;
00248           ev.getByLabel( theMatchedHits, matchedHits);
00249           std::pair<DetId,DetIdTECSameDiskComparator> getter = accessor.stripTECDisk(theSide,theIdLayer);
00250           SiStripMatchedRecHit2DCollection::Range range = matchedHits->equal_range(getter.first, getter.second);
00251           for (SiStripMatchedRecHit2DCollection::const_iterator it = range.first; it != range.second; ++it) {
00252               int ring = TECDetId( it->detId() ).ring();  if (!ringRange(ring)) continue;
00253               for (SiStripMatchedRecHit2DCollection::DetSet::const_iterator hit = it->begin(), end = it->end(); hit != end; ++hit) {
00254                   result.push_back(  sl.hitBuilder()->build(hit) );
00255               }
00256           }
00257           if (skipClusters) cleanedOfClusters(ev,result,true);
00258       }
00259       if (hasRPhiHits) {
00260           edm::Handle<SiStripRecHit2DCollection> rphiHits;
00261           ev.getByLabel( theRPhiHits, rphiHits);
00262           std::pair<DetId,DetIdTECSameDiskComparator> getter = accessor.stripTECDisk(theSide,theIdLayer);
00263           SiStripRecHit2DCollection::Range range = rphiHits->equal_range(getter.first, getter.second);
00264           for (SiStripRecHit2DCollection::const_iterator it = range.first; it != range.second; ++it) {
00265               int ring = TECDetId( it->detId() ).ring();  if (!ringRange(ring)) continue;
00266               if ((SiStripDetId(it->detId()).partnerDetId() != 0) && hasSimpleRphiHitsCleaner) continue;  // this is a brutal "cleaning". Add something smarter in the future
00267               for (SiStripRecHit2DCollection::DetSet::const_iterator hit = it->begin(), end = it->end(); hit != end; ++hit) {
00268                   result.push_back( sl.hitBuilder()->build(hit) );
00269               }
00270           }
00271           if (skipClusters) cleanedOfClusters(ev,result,false);
00272 
00273       }
00274       if (hasStereoHits) {
00275           edm::Handle<SiStripRecHit2DCollection> stereoHits;
00276           ev.getByLabel( theStereoHits, stereoHits);
00277           std::pair<DetId,DetIdTECSameDiskComparator> getter = accessor.stripTECDisk(theSide,theIdLayer);
00278           SiStripRecHit2DCollection::Range range = stereoHits->equal_range(getter.first, getter.second);
00279           for (SiStripRecHit2DCollection::const_iterator it = range.first; it != range.second; ++it) {
00280               int ring = TECDetId( it->detId() ).ring();  if (!ringRange(ring)) continue;
00281               for (SiStripRecHit2DCollection::DetSet::const_iterator hit = it->begin(), end = it->end(); hit != end; ++hit) {
00282                   result.push_back( sl.hitBuilder()->build(hit) );
00283               }
00284           }
00285           if (skipClusters) cleanedOfClusters(ev,result,false);
00286       }
00287   }
00288   LogDebug("HitExtractorSTRP")<<" giving: "<<result.size()<<" out";
00289   return result;
00290 }
00291 
00292