CMS 3D CMS Logo

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