CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch13/src/RecoLocalTracker/SiStripRecHitConverter/src/SiStripRecHitConverterAlgorithm.cc

Go to the documentation of this file.
00001 #include "RecoLocalTracker/SiStripRecHitConverter/interface/SiStripRecHitConverterAlgorithm.h"
00002 #include "RecoLocalTracker/Records/interface/TkStripCPERecord.h"
00003 #include "CalibTracker/Records/interface/SiStripQualityRcd.h"
00004 
00005 #include "Geometry/TrackerGeometryBuilder/interface/GluedGeomDet.h"
00006 #include "Geometry/TrackerGeometryBuilder/interface/StripGeomDetUnit.h"
00007 
00008 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
00009 #include "DataFormats/Common/interface/RefGetter.h"
00010 #include "DataFormats/Common/interface/Ref.h"
00011 
00012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00013 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00014 #include "FWCore/Utilities/interface/typelookup.h"
00015 TYPELOOKUP_DATA_REG(SiStripRecHitMatcher);
00016 
00017 
00018 SiStripRecHitConverterAlgorithm::SiStripRecHitConverterAlgorithm(const edm::ParameterSet& conf) : 
00019   useQuality(conf.getParameter<bool>("useSiStripQuality")),
00020   maskBad128StripBlocks( conf.existsAs<bool>("MaskBadAPVFibers") && conf.getParameter<bool>("MaskBadAPVFibers")),
00021   tracker_cache_id(0),
00022   cpe_cache_id(0),
00023   quality_cache_id(0),
00024   cpeTag(conf.getParameter<edm::ESInputTag>("StripCPE")),
00025   matcherTag(conf.getParameter<edm::ESInputTag>("Matcher")),
00026   qualityTag(conf.getParameter<edm::ESInputTag>("siStripQualityLabel"))
00027 {}
00028 
00029 void SiStripRecHitConverterAlgorithm::
00030 initialize(const edm::EventSetup& es) 
00031 {
00032   uint32_t tk_cache_id = es.get<TrackerDigiGeometryRecord>().cacheIdentifier();
00033   uint32_t c_cache_id = es.get<TkStripCPERecord>().cacheIdentifier();
00034   uint32_t q_cache_id = es.get<SiStripQualityRcd>().cacheIdentifier();
00035 
00036   if(tk_cache_id != tracker_cache_id) {
00037     es.get<TrackerDigiGeometryRecord>().get(tracker);
00038     tracker_cache_id = tk_cache_id;
00039   }
00040   if(c_cache_id != cpe_cache_id) {
00041     es.get<TkStripCPERecord>().get(matcherTag, matcher);
00042     es.get<TkStripCPERecord>().get(cpeTag, parameterestimator); 
00043     cpe_cache_id = c_cache_id;
00044   }
00045   if( useQuality && q_cache_id!=quality_cache_id) {
00046     es.get<SiStripQualityRcd>().get(qualityTag, quality);
00047     quality_cache_id = q_cache_id;
00048   }
00049 }
00050 
00051 void SiStripRecHitConverterAlgorithm::
00052 run(edm::Handle<edmNew::DetSetVector<SiStripCluster> > input, products& output) 
00053 { run(input, output, LocalVector(0.,0.,0.)); }
00054 
00055 void SiStripRecHitConverterAlgorithm::
00056 run(edm::Handle<edmNew::DetSetVector<SiStripCluster> > inputhandle, products& output, LocalVector trackdirection)
00057 {
00058   for (edmNew::DetSetVector<SiStripCluster>::const_iterator 
00059          DS = inputhandle->begin(); DS != inputhandle->end(); DS++ ) {     if(!useModule(DS->id())) continue;
00060 
00061     Collector collector = StripSubdetector(DS->id()).stereo()  
00062       ? Collector(*output.stereo, DS->id()) 
00063       : Collector(*output.rphi,   DS->id());
00064     bool bad128StripBlocks[6]; fillBad128StripBlocks( DS->id(), bad128StripBlocks);
00065     
00066 
00067     for(edmNew::DetSet<SiStripCluster>::const_iterator 
00068           cluster = DS->begin();  cluster != DS->end(); ++cluster ) {     if(isMasked(*cluster,bad128StripBlocks)) continue;
00069 
00070       StripClusterParameterEstimator::LocalValues parameters =  parameterestimator->localParameters(*cluster,*(tracker->idToDetUnit(DS->id())));
00071       collector.push_back(SiStripRecHit2D( parameters.first, parameters.second, DS->id(), edmNew::makeRefTo(inputhandle,cluster) ));
00072     }
00073 
00074     if (collector.empty()) collector.abort();
00075   }
00076   match(output,trackdirection);
00077 }
00078 
00079 void SiStripRecHitConverterAlgorithm::
00080 run(edm::Handle<edm::RefGetter<SiStripCluster> >  refGetterhandle, 
00081     edm::Handle<edm::LazyGetter<SiStripCluster> >  lazyGetterhandle, 
00082     products& output)
00083 {
00084   bool bad128StripBlocks[6];
00085   bool goodDet=true;
00086   DetId lastId(0);
00087   std::auto_ptr<Collector> collector(new Collector(*output.stereo, lastId));
00088   
00089   edm::RefGetter<SiStripCluster>::const_iterator iregion = refGetterhandle->begin();
00090   for(;iregion!=refGetterhandle->end();++iregion) {
00091     const edm::RegionIndex<SiStripCluster>& region = *iregion;
00092     const uint32_t start = region.start();
00093     const uint32_t finish = region.finish();
00094     for (uint32_t i = start; i < finish; i++) {
00095       edm::RegionIndex<SiStripCluster>::const_iterator icluster = region.begin()+(i-start);
00096       
00097       DetId detId(icluster->geographicalId());
00098       if(detId != lastId) {
00099         if(collector->empty()) collector->abort();
00100         lastId = detId;
00101         goodDet = useModule(detId);
00102         if(goodDet) {
00103           fillBad128StripBlocks(detId, bad128StripBlocks);
00104           collector = StripSubdetector(detId).stereo()  
00105             ? std::auto_ptr<Collector>(new Collector(*output.stereo, detId)) 
00106             : std::auto_ptr<Collector>(new Collector(*output.rphi,   detId));
00107         }
00108       }
00109       if( !goodDet || isMasked(*icluster, bad128StripBlocks)) continue;
00110 
00111       StripClusterParameterEstimator::LocalValues parameters = parameterestimator->localParameters(*icluster,*(tracker->idToDetUnit(detId)));
00112       collector->push_back(SiStripRecHit2D( parameters.first, parameters.second, detId, makeRefToLazyGetter(lazyGetterhandle,i) ));      
00113     }
00114     if(collector->empty()) collector->abort();
00115   }
00116   match(output,LocalVector(0.,0.,0.));
00117 }
00118 
00119 
00120 namespace {
00121 
00122   struct CollectorHelper {
00123     bool regional;
00124     size_t nmatch;
00125     
00126     
00127     typedef edm::OwnVector<SiStripMatchedRecHit2D> CollectorMatched;
00128     typedef SiStripMatchedRecHit2DCollection::FastFiller Collector;
00129     
00130     Collector & m_collector;
00131     CollectorMatched & m_collectorMatched;
00132     SiStripRecHit2DCollection::FastFiller & m_fillerRphiUnm;
00133     std::vector<SiStripRecHit2D::ClusterRef::key_type>         & m_matchedSteroClusters;
00134     std::vector<SiStripRecHit2D::ClusterRegionalRef::key_type> & m_matchedSteroClustersRegional;
00135     
00136     static inline SiStripRecHit2D const & stereoHit(edmNew::DetSet<SiStripRecHit2D>::const_iterator iter) {
00137       return *iter;
00138     }
00139     
00140     static inline SiStripRecHit2D const & monoHit(edmNew::DetSet<SiStripRecHit2D>::const_iterator iter) {
00141       return *iter;
00142     }
00143     
00144     struct Add {
00145       Add(CollectorHelper& ih) : h(ih){}
00146       CollectorHelper& h;
00147       void operator()(SiStripMatchedRecHit2D const & rh) { h.m_collectorMatched.push_back(rh);}
00148     };
00149 
00150     CollectorHelper & collector() {
00151       return *this;
00152     }
00153 
00154    void operator()(SiStripMatchedRecHit2D const & rh) {m_collectorMatched.push_back(rh);}
00155 
00156 
00157     CollectorHelper(
00158                     Collector & i_collector,
00159                     CollectorMatched & i_collectorMatched,
00160                     SiStripRecHit2DCollection::FastFiller & i_fillerRphiUnm,
00161                     std::vector<SiStripRecHit2D::ClusterRef::key_type>         & i_matchedSteroClusters,
00162                     std::vector<SiStripRecHit2D::ClusterRegionalRef::key_type> & i_matchedSteroClustersRegional
00163                     ) : regional(false), nmatch(0), 
00164                         m_collector(i_collector),
00165                         m_collectorMatched(i_collectorMatched),
00166                         m_fillerRphiUnm(i_fillerRphiUnm),
00167                         m_matchedSteroClusters(i_matchedSteroClusters),
00168                         m_matchedSteroClustersRegional(i_matchedSteroClustersRegional)
00169     {}
00170     
00171     void closure(edmNew::DetSet<SiStripRecHit2D>::const_iterator it) {
00172       if (!m_collectorMatched.empty()){
00173         nmatch+=m_collectorMatched.size();
00174         for (edm::OwnVector<SiStripMatchedRecHit2D>::const_iterator itm = m_collectorMatched.begin(),
00175                edm = m_collectorMatched.end();
00176              itm != edm; 
00177              ++itm) {
00178           m_collector.push_back(*itm);
00179           // mark the stereo hit cluster as used, so that the hit won't go in the unmatched stereo ones
00180           if (itm->stereoHit()->cluster().isNonnull()) {
00181             m_matchedSteroClusters.push_back(itm->stereoHit()->cluster().key()); 
00182           } else {
00183             m_matchedSteroClustersRegional.push_back(itm->stereoHit()->cluster_regional().key()); 
00184             regional = true;
00185           }
00186         }
00187         m_collectorMatched.clear();
00188       } else {
00189         // store a copy of this rphi hit as an unmatched rphi hit
00190         m_fillerRphiUnm.push_back(*it);
00191       }
00192     }
00193   };
00194 }
00195 
00196 
00197 void SiStripRecHitConverterAlgorithm::
00198 match(products& output, LocalVector trackdirection) const 
00199 {
00200   int nmatch=0;
00201   edm::OwnVector<SiStripMatchedRecHit2D> collectorMatched; // gp/FIXME: avoid this
00202   
00203   // Remember the ends of the collections, as we will use them a lot
00204   SiStripRecHit2DCollection::const_iterator edStereoDet = output.stereo->end();
00205   SiStripRecHit2DCollection::const_iterator edRPhiDet   = output.rphi->end();
00206   
00207   // two work vectors for bookeeping clusters used by the stereo part of the matched hits
00208   std::vector<SiStripRecHit2D::ClusterRef::key_type>         matchedSteroClusters;
00209   std::vector<SiStripRecHit2D::ClusterRegionalRef::key_type> matchedSteroClustersRegional;
00210   
00211   for (SiStripRecHit2DCollection::const_iterator itRPhiDet = output.rphi->begin(); itRPhiDet != edRPhiDet; ++itRPhiDet) {
00212     edmNew::DetSet<SiStripRecHit2D> rphiHits = *itRPhiDet;
00213     StripSubdetector specDetId(rphiHits.detId());
00214     uint32_t partnerId = specDetId.partnerDetId();
00215     
00216     // if not part of a glued pair
00217     if (partnerId == 0) { 
00218       // I must copy these as unmatched 
00219       if (!rphiHits.empty()) {
00220         SiStripRecHit2DCollection::FastFiller filler(*output.rphiUnmatched, rphiHits.detId());
00221         filler.resize(rphiHits.size());
00222         std::copy(rphiHits.begin(), rphiHits.end(), filler.begin());
00223       }
00224       continue;
00225     }
00226     
00227     SiStripRecHit2DCollection::const_iterator itStereoDet = output.stereo->find(partnerId);
00228     
00229     // if the partner is not found (which probably can happen if it's empty)
00230     if (itStereoDet == edStereoDet) {
00231       // I must copy these as unmatched 
00232       if (!rphiHits.empty()) {
00233         SiStripRecHit2DCollection::FastFiller filler(*output.rphiUnmatched, rphiHits.detId());
00234         filler.resize(rphiHits.size());
00235         std::copy(rphiHits.begin(), rphiHits.end(), filler.begin());
00236       }
00237       continue;
00238     }
00239     
00240     edmNew::DetSet<SiStripRecHit2D> stereoHits = *itStereoDet;
00241     
00242      
00243     // Get ready for making glued hits
00244     const GluedGeomDet* gluedDet = (const GluedGeomDet*)tracker->idToDet(DetId(specDetId.glued()));
00245     typedef SiStripMatchedRecHit2DCollection::FastFiller Collector;
00246     Collector collector(*output.matched, specDetId.glued());
00247     
00248     // Prepare also the list for unmatched rphi hits
00249     SiStripRecHit2DCollection::FastFiller fillerRphiUnm(*output.rphiUnmatched, rphiHits.detId());
00250     
00251     // a list of clusters used by the matched part of the stereo hits in this detector
00252     matchedSteroClusters.clear();          // at the beginning, empty
00253     matchedSteroClustersRegional.clear();  // I need two because the refs can be different
00254     bool regional = false;                 // I also want to remember if they come from standard or HLT reco
00255 
00256 #ifdef DOUBLE_MATCH
00257     CollectorHelper chelper(collector, collectorMatched,
00258                     fillerRphiUnm,
00259                     matchedSteroClusters,matchedSteroClustersRegional
00260                     );
00261     matcher->doubleMatch(rphiHits.begin(), rphiHits.end(), 
00262                          stereoHits.begin(),stereoHits.end(),gluedDet,trackdirection,chelper);
00263     regional = chelper.regional;
00264     nmatch+=chelper.nmatch;
00265 #else 
00266    // Make simple collection of this (gp:FIXME: why do we need it?)
00267     SiStripRecHitMatcher::SimpleHitCollection stereoSimpleHits;
00268     // gp:FIXME: use std::transform 
00269     stereoSimpleHits.reserve(stereoHits.size());
00270     for (edmNew::DetSet<SiStripRecHit2D>::const_iterator it = stereoHits.begin(), ed = stereoHits.end(); it != ed; ++it) {
00271       stereoSimpleHits.push_back(&*it);
00272     }
00273 
00274     for (edmNew::DetSet<SiStripRecHit2D>::const_iterator it = rphiHits.begin(), ed = rphiHits.end(); it != ed; ++it) {
00275       matcher->match(&(*it),stereoSimpleHits.begin(),stereoSimpleHits.end(),collectorMatched,gluedDet,trackdirection);
00276       if (collectorMatched.size()>0){
00277         nmatch+=collectorMatched.size();
00278         for (edm::OwnVector<SiStripMatchedRecHit2D>::const_iterator itm = collectorMatched.begin(),
00279                edm = collectorMatched.end();
00280              itm != edm; 
00281              ++itm) {
00282           collector.push_back(*itm);
00283           // mark the stereo hit cluster as used, so that the hit won't go in the unmatched stereo ones
00284           if (itm->stereoHit()->cluster().isNonnull()) {
00285             matchedSteroClusters.push_back(itm->stereoHit()->cluster().key()); 
00286           } else {
00287             matchedSteroClustersRegional.push_back(itm->stereoHit()->cluster_regional().key()); 
00288             regional = true;
00289           }
00290         }
00291         collectorMatched.clear();
00292       } else {
00293         // store a copy of this rphi hit as an unmatched rphi hit
00294         fillerRphiUnm.push_back(*it);
00295       }
00296     }
00297     
00298 #endif
00299 
00300 
00301     // discard matched hits if the collection is empty
00302     if (collector.empty()) collector.abort();
00303     
00304     // discard unmatched rphi hits if there are none
00305     if (fillerRphiUnm.empty()) fillerRphiUnm.abort();
00306     
00307     // now look for unmatched stereo hits    
00308     SiStripRecHit2DCollection::FastFiller fillerStereoUnm(*output.stereoUnmatched, stereoHits.detId());
00309     if (!regional) {
00310       std::sort(matchedSteroClusters.begin(), matchedSteroClusters.end());
00311       for (edmNew::DetSet<SiStripRecHit2D>::const_iterator it = stereoHits.begin(), ed = stereoHits.end(); it != ed; ++it) {
00312         if (!std::binary_search(matchedSteroClusters.begin(), matchedSteroClusters.end(), it->cluster().key())) {
00313           fillerStereoUnm.push_back(*it);
00314         }
00315       }
00316     } else {
00317       std::sort(matchedSteroClustersRegional.begin(), matchedSteroClustersRegional.end());
00318       for (edmNew::DetSet<SiStripRecHit2D>::const_iterator it = stereoHits.begin(), ed = stereoHits.end(); it != ed; ++it) {
00319         if (!std::binary_search(matchedSteroClustersRegional.begin(), matchedSteroClustersRegional.end(), it->cluster_regional().key())) {
00320           fillerStereoUnm.push_back(*it);
00321         }
00322       }
00323     }
00324     if (fillerStereoUnm.empty()) fillerStereoUnm.abort(); 
00325     
00326     
00327   }
00328   
00329   for (SiStripRecHit2DCollection::const_iterator itStereoDet = output.stereo->begin(); itStereoDet != edStereoDet; ++itStereoDet) {
00330     edmNew::DetSet<SiStripRecHit2D> stereoHits = *itStereoDet;
00331     StripSubdetector specDetId(stereoHits.detId());
00332     uint32_t partnerId = specDetId.partnerDetId();
00333     if (partnerId == 0) continue;
00334     SiStripRecHit2DCollection::const_iterator itRPhiDet = output.rphi->find(partnerId);
00335     if (itRPhiDet == edRPhiDet) {
00336       if (!stereoHits.empty()) {
00337         SiStripRecHit2DCollection::FastFiller filler(*output.stereoUnmatched, stereoHits.detId());
00338         filler.resize(stereoHits.size());
00339         std::copy(stereoHits.begin(), stereoHits.end(), filler.begin());
00340       }
00341     }
00342   }
00343   
00344   edm::LogInfo("SiStripRecHitConverter") 
00345     << "found\n"         
00346     << nmatch 
00347     << "  matched RecHits\n";
00348 }
00349 
00350 void SiStripRecHitConverterAlgorithm::
00351 fillBad128StripBlocks(const uint32_t detid, bool bad128StripBlocks[6] ) const 
00352 {
00353   if(maskBad128StripBlocks) {
00354     short badApvs   = quality->getBadApvs(detid);
00355     short badFibers = quality->getBadFibers(detid);
00356     for (int j = 0; j < 6; j++) {
00357       bad128StripBlocks[j] = (badApvs & (1 << j));
00358     }
00359     for (int j = 0; j < 3; j++) {
00360       if (badFibers & (1 << j)) {
00361         bad128StripBlocks[2*j+0] = true;
00362         bad128StripBlocks[2*j+1] = true;
00363       }
00364     }
00365   }
00366 }
00367 
00368 inline
00369 bool SiStripRecHitConverterAlgorithm::
00370 isMasked(const SiStripCluster &cluster, bool bad128StripBlocks[6]) const 
00371 {
00372   if(maskBad128StripBlocks) {
00373     if ( bad128StripBlocks[cluster.firstStrip() >> 7] ) {
00374       if ( bad128StripBlocks[(cluster.firstStrip()+cluster.amplitudes().size())  >> 7] ||
00375            bad128StripBlocks[static_cast<int32_t>(cluster.barycenter()-0.499999) >> 7] ) {
00376         return true;
00377       }
00378     } else {
00379       if ( bad128StripBlocks[(cluster.firstStrip()+cluster.amplitudes().size())  >> 7] &&
00380            bad128StripBlocks[static_cast<int32_t>(cluster.barycenter()-0.499999) >> 7] ) {
00381         return true;
00382       }
00383     }
00384   }
00385   return false;
00386 }
00387 
00388 inline
00389 bool SiStripRecHitConverterAlgorithm::
00390 useModule(const uint32_t id) const
00391 {
00392   const StripGeomDetUnit * stripdet=(const StripGeomDetUnit*)tracker->idToDetUnit(id);
00393   if(stripdet==0) edm::LogWarning("SiStripRecHitConverter") << "Detid=" << id << " not found";
00394   return stripdet!=0 && (!useQuality || quality->IsModuleUsable(id));
00395 }