CMS 3D CMS Logo

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