CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_13_patch3/src/RecoTracker/MeasurementDet/plugins/TkStripMeasurementDet.cc

Go to the documentation of this file.
00001 #include "TkStripMeasurementDet.h"
00002 #include "TrackingTools/TransientTrackingRecHit/interface/InvalidTransientRecHit.h"
00003 #include "Geometry/CommonTopologies/interface/StripTopology.h"
00004 #include "TrackingTools/MeasurementDet/interface/MeasurementDetException.h"
00005 #include "TrackingTools/PatternTools/interface/TrajectoryMeasurement.h"
00006 #include "StripClusterAboveU.h"
00007 #include "RecoTracker/TransientTrackingRecHit/interface/TSiStripRecHit2DLocalPos.h"
00008 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2D.h"
00009 #include "TrackingTools/PatternTools/interface/MeasurementEstimator.h"
00010 #include "TrackingTools/PatternTools/interface/TrajMeasLessEstim.h"
00011 
00012 #include <typeinfo>
00013 // #include "TrackingTools/KalmanUpdators/interface/Chi2MeasurementEstimator.h"
00014 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00015 #include "DataFormats/SiStripDetId/interface/SiStripDetId.h"
00016 
00017 TkStripMeasurementDet::TkStripMeasurementDet( const GeomDet* gdet,
00018                                               StMeasurementDetSet & dets
00019                                               ) : 
00020   MeasurementDet (gdet),theDets_(&dets), index_(-1)
00021   {
00022     if (dynamic_cast<const StripGeomDetUnit*>(gdet) == 0) {
00023       throw MeasurementDetException( "TkStripMeasurementDet constructed with a GeomDet which is not a StripGeomDetUnit");
00024     }
00025   }
00026 
00027 
00028 std::vector<TrajectoryMeasurement> 
00029 TkStripMeasurementDet::
00030 fastMeasurements( const TrajectoryStateOnSurface& stateOnThisDet, 
00031                   const TrajectoryStateOnSurface& startingState, 
00032                   const Propagator&, 
00033                   const MeasurementEstimator& est) const
00034 { 
00035   std::vector<TrajectoryMeasurement> result;
00036 
00037   if (!isActive()) {
00038     LogDebug("TkStripMeasurementDet")<<" found an inactive module "<<rawId();
00039     result.push_back( TrajectoryMeasurement( stateOnThisDet, 
00040                 InvalidTransientRecHit::build(&geomDet(), TrackingRecHit::inactive), 
00041                 0.F));
00042     return result;
00043   }
00044  
00045   float utraj =  specificGeomDet().specificTopology().measurementPosition( stateOnThisDet.localPosition()).x();
00046   float uerr;
00047   //  if (theClusterRange.first == theClusterRange.second) { // empty
00048   if (isEmpty()){
00049     LogDebug("TkStripMeasurementDet") << " DetID " << rawId() << " empty ";
00050     if (stateOnThisDet.hasError()){
00051     uerr= sqrt(specificGeomDet().specificTopology().measurementError(stateOnThisDet.localPosition(),stateOnThisDet.localError().positionError()).uu());
00052      if (testStrips(utraj,uerr)) {
00053         result.push_back( TrajectoryMeasurement( stateOnThisDet, InvalidTransientRecHit::build(&fastGeomDet(), TrackingRecHit::missing), 0.F));
00054      } else { 
00055         result.push_back( TrajectoryMeasurement( stateOnThisDet, InvalidTransientRecHit::build(&fastGeomDet(), TrackingRecHit::inactive), 0.F));
00056      }
00057     }else{
00058       result.push_back( TrajectoryMeasurement( stateOnThisDet, InvalidTransientRecHit::build(&fastGeomDet(), TrackingRecHit::missing), 0.F));
00059     }
00060     return result;
00061   }
00062   
00063   if(!isRegional()){//old implemetation with DetSet
00064     auto rightCluster = 
00065       std::find_if( detSet().begin(), detSet().end(), StripClusterAboveU( utraj)); //FIXME
00066     
00067     if ( rightCluster != detSet().begin()) {
00068       // there are hits on the left of the utraj
00069       new_const_iterator leftCluster = rightCluster;
00070       while ( --leftCluster >=  detSet().begin()) {
00071         if (isMasked(*leftCluster)) continue;
00072         SiStripClusterRef clusterref = edmNew::makeRefTo( handle(), leftCluster ); 
00073         if (accept(clusterref)){
00074           RecHitContainer recHits = buildRecHits(clusterref,stateOnThisDet); 
00075           bool isCompatible(false);
00076           for(RecHitContainer::const_iterator recHit=recHits.begin();recHit!=recHits.end();++recHit){     
00077             std::pair<bool,double> diffEst = est.estimate(stateOnThisDet, **recHit);
00078             if ( diffEst.first ) {
00079               result.push_back( TrajectoryMeasurement( stateOnThisDet, *recHit, 
00080                                                        diffEst.second));
00081               isCompatible = true;
00082             }
00083           }
00084           if(!isCompatible) break; // exit loop on first incompatible hit
00085         }
00086         else LogDebug("TkStripMeasurementDet")<<"skipping this str from last iteration on"<<rawId()<<" key: "<<clusterref.key();
00087       }
00088     }
00089     
00090 
00091     for ( ; rightCluster != detSet().end(); rightCluster++) {
00092       if (isMasked(*rightCluster)) continue;
00093       SiStripClusterRef clusterref = edmNew::makeRefTo( handle(), rightCluster ); 
00094       if (accept(clusterref)){
00095         RecHitContainer recHits = buildRecHits(clusterref,stateOnThisDet); 
00096         bool isCompatible(false);
00097         for(RecHitContainer::const_iterator recHit=recHits.begin();recHit!=recHits.end();++recHit){       
00098           std::pair<bool,double> diffEst = est.estimate(stateOnThisDet, **recHit);
00099           if ( diffEst.first ) {
00100             result.push_back( TrajectoryMeasurement( stateOnThisDet, *recHit, 
00101                                                      diffEst.second));
00102             isCompatible = true;
00103           }
00104         }
00105         if(!isCompatible) break; // exit loop on first incompatible hit
00106       }
00107       else LogDebug("TkStripMeasurementDet")<<"skipping this str from last iteration on" << rawId()<<" key: "<<clusterref.key();
00108     }
00109   }// end block with DetSet
00110   else{
00111     LogDebug("TkStripMeasurementDet")<<" finding left/ right";
00112     result.reserve(size());
00113     unsigned int rightCluster = beginClusterI();
00114     for (; rightCluster!= endClusterI();++rightCluster){
00115       SiStripRegionalClusterRef clusterref = edm::makeRefToLazyGetter(regionalHandle(),rightCluster);
00116       if (clusterref->barycenter() > utraj) break;
00117     }
00118     
00119     unsigned int leftCluster = 1;
00120     for (unsigned int iReadBackWard=1; iReadBackWard<=(rightCluster-beginClusterI()) ; ++iReadBackWard){
00121       leftCluster=rightCluster-iReadBackWard;
00122       SiStripRegionalClusterRef clusterref = edm::makeRefToLazyGetter(regionalHandle(),leftCluster);
00123       if (isMasked(*clusterref)) continue;
00124       if (accept(clusterref)){
00125         RecHitContainer recHits = buildRecHits(clusterref,stateOnThisDet); 
00126         bool isCompatible(false);
00127         for(RecHitContainer::const_iterator recHit=recHits.begin();recHit!=recHits.end();++recHit){       
00128           std::pair<bool,double> diffEst = est.estimate(stateOnThisDet, **recHit);
00129           if ( diffEst.first ) {
00130             result.push_back( TrajectoryMeasurement( stateOnThisDet, *recHit, 
00131                                                      diffEst.second));
00132             isCompatible = true;
00133           }
00134         }
00135         if(!isCompatible) break; // exit loop on first incompatible hit
00136       }
00137       else LogDebug("TkStripMeasurementDet")<<"skipping this reg str from last iteration on"<<rawId()<<" key: "<<clusterref.key();
00138     }
00139     
00140     
00141     for ( ; rightCluster != endClusterI(); ++rightCluster) {
00142       SiStripRegionalClusterRef clusterref = edm::makeRefToLazyGetter(regionalHandle(),rightCluster);
00143       if (isMasked(*clusterref)) continue;
00144       if (accept(clusterref)){
00145         RecHitContainer recHits = buildRecHits(clusterref,stateOnThisDet); 
00146         bool isCompatible(false);
00147         for(RecHitContainer::const_iterator recHit=recHits.begin();recHit!=recHits.end();++recHit){       
00148           std::pair<bool,double> diffEst = est.estimate(stateOnThisDet, **recHit);
00149           if ( diffEst.first ) {
00150             result.push_back( TrajectoryMeasurement( stateOnThisDet, *recHit, 
00151                                                      diffEst.second));
00152             isCompatible = true;
00153           }
00154         }
00155         if(!isCompatible) break; // exit loop on first incompatible hit
00156       }
00157       else LogDebug("TkStripMeasurementDet")<<"skipping this reg str from last iteration on"<<rawId()<<" key: "<<clusterref.key();
00158     }
00159   }
00160   
00161 
00162   if ( result.empty()) {
00163     // create a TrajectoryMeasurement with an invalid RecHit and zero estimate
00164     if (stateOnThisDet.hasError()){
00165     uerr= sqrt(specificGeomDet().specificTopology().measurementError(stateOnThisDet.localPosition(),stateOnThisDet.localError().positionError()).uu());
00166      if (testStrips(utraj,uerr)) {
00167        //LogDebug("TkStripMeasurementDet") << " DetID " << id_ << " empty after search, but active ";
00168        result.push_back( TrajectoryMeasurement( stateOnThisDet, InvalidTransientRecHit::build(&fastGeomDet(), TrackingRecHit::missing), 0.F));
00169      } else { 
00170        //LogDebug("TkStripMeasurementDet") << " DetID " << id_ << " empty after search, and inactive ";
00171        result.push_back( TrajectoryMeasurement( stateOnThisDet, InvalidTransientRecHit::build(&fastGeomDet(), TrackingRecHit::inactive), 0.F));
00172      }
00173     }else{
00174       result.push_back( TrajectoryMeasurement( stateOnThisDet, InvalidTransientRecHit::build(&fastGeomDet(), TrackingRecHit::missing), 0.F));
00175     }
00176   }
00177   else {
00178     //LogDebug("TkStripMeasurementDet") << " DetID " << id_ << " full: " << (result.size()) << " compatible hits";
00179     // sort results according to estimator value
00180     if ( result.size() > 1) {
00181       std::sort( result.begin(), result.end(), TrajMeasLessEstim());
00182     }
00183   }
00184   return result;
00185 }
00186 
00187 
00188 TransientTrackingRecHit::RecHitPointer
00189 TkStripMeasurementDet::buildRecHit( const SiStripClusterRef& cluster,
00190                                     const TrajectoryStateOnSurface& ltp) const
00191 {
00192   const GeomDetUnit& gdu( specificGeomDet());
00193   LocalValues lv = cpe()->localParameters( *cluster, gdu, ltp);
00194   return TSiStripRecHit2DLocalPos::build( lv.first, lv.second, &fastGeomDet(), cluster, cpe());
00195 }
00196 
00197 TransientTrackingRecHit::RecHitPointer
00198 TkStripMeasurementDet::buildRecHit( const SiStripRegionalClusterRef& cluster,
00199                                     const TrajectoryStateOnSurface& ltp) const
00200 {
00201   const GeomDetUnit& gdu( specificGeomDet());
00202   LocalValues lv = cpe()->localParameters( *cluster, gdu, ltp);
00203   return TSiStripRecHit2DLocalPos::build( lv.first, lv.second, &fastGeomDet(), cluster, cpe());
00204 }
00205 
00206 
00207 
00208 TkStripMeasurementDet::RecHitContainer 
00209 TkStripMeasurementDet::buildRecHits( const SiStripClusterRef& cluster,
00210                                      const TrajectoryStateOnSurface& ltp) const
00211 {
00212   const GeomDetUnit& gdu( specificGeomDet());
00213   VLocalValues vlv = cpe()->localParametersV( *cluster, gdu, ltp);
00214   RecHitContainer res;
00215   for(VLocalValues::const_iterator it=vlv.begin();it!=vlv.end();++it){
00216     res.push_back(TSiStripRecHit2DLocalPos::build( it->first, it->second, &fastGeomDet(), cluster, cpe()));
00217   }
00218   return res; 
00219 }
00220 
00221 
00222 TkStripMeasurementDet::RecHitContainer 
00223 TkStripMeasurementDet::buildRecHits( const SiStripRegionalClusterRef& cluster,
00224                                      const TrajectoryStateOnSurface& ltp) const
00225 {
00226   const GeomDetUnit& gdu( specificGeomDet());
00227   VLocalValues vlv = cpe()->localParametersV( *cluster, gdu, ltp);
00228   RecHitContainer res;
00229   for(VLocalValues::const_iterator it=vlv.begin();it!=vlv.end();++it){
00230     res.push_back(TSiStripRecHit2DLocalPos::build( it->first, it->second, &fastGeomDet(), cluster, cpe()));
00231   }
00232   return res; 
00233 }
00234 
00235 
00236 
00237 TkStripMeasurementDet::RecHitContainer 
00238 TkStripMeasurementDet::recHits( const TrajectoryStateOnSurface& ts) const
00239 {
00240   RecHitContainer result;
00241   if (isEmpty() == true) return result;
00242   if (isActive() == false) return result; // GIO
00243 
00244   if(!isRegional()){//old implemetation with DetSet
00245     result.reserve(detSet().size());
00246     for ( new_const_iterator ci = detSet().begin(); ci != detSet().end(); ++ ci ) {
00247       if (isMasked(*ci)) continue;
00248       // for ( ClusterIterator ci=theClusterRange.first; ci != theClusterRange.second; ci++) {
00249       SiStripClusterRef  cluster = edmNew::makeRefTo( handle(), ci ); 
00250       if (accept(cluster))
00251         result.push_back( buildRecHit( cluster, ts));
00252       else LogDebug("TkStripMeasurementDet")<<"skipping this str from last iteration on"<<rawId()<<" key: "<<cluster.key();
00253     }
00254   }else{
00255     result.reserve(size());
00256     for (unsigned int ci = beginClusterI() ; ci!= endClusterI();++ci){
00257       SiStripRegionalClusterRef clusterRef = edm::makeRefToLazyGetter(regionalHandle(),ci);
00258       if (isMasked(*clusterRef)) continue;
00259       if (accept(clusterRef))
00260         result.push_back( buildRecHit( clusterRef, ts));
00261       else LogDebug("TkStripMeasurementDet")<<"skipping this reg str from last iteration on"<<rawId()<<" key: "<<clusterRef.key();
00262       }
00263   }
00264   return result;
00265 
00266 }
00267 
00268 template<class ClusterRefT>
00269 void
00270 TkStripMeasurementDet::buildSimpleRecHit( const ClusterRefT& cluster,
00271                                           const TrajectoryStateOnSurface& ltp,
00272                                           std::vector<SiStripRecHit2D>& res ) const
00273 {
00274   const GeomDetUnit& gdu( specificGeomDet());
00275   VLocalValues vlv = cpe()->localParametersV( *cluster, gdu, ltp);
00276   for(VLocalValues::const_iterator it=vlv.begin();it!=vlv.end();++it){
00277     res.push_back(SiStripRecHit2D( it->first, it->second, rawId(), cluster));
00278   }
00279 }
00280 
00281 
00282 void 
00283 TkStripMeasurementDet::simpleRecHits( const TrajectoryStateOnSurface& ts, std::vector<SiStripRecHit2D> &result) const
00284 {
00285   if (isEmpty() || !isActive()) return;
00286 
00287   if(!isRegional()){//old implemetation with DetSet
00288     result.reserve(detSet().size());
00289     for ( new_const_iterator ci = detSet().begin(); ci != detSet().end(); ++ ci ) {
00290       if (isMasked(*ci)) continue;
00291       // for ( ClusterIterator ci=theClusterRange.first; ci != theClusterRange.second; ci++) {
00292       SiStripClusterRef  cluster = edmNew::makeRefTo( handle(), ci ); 
00293       if (accept(cluster))
00294         buildSimpleRecHit( cluster, ts,result);
00295       else LogDebug("TkStripMeasurementDet")<<"skipping this str from last iteration on"<<rawId()<<" key: "<<cluster.key();
00296     }
00297   }else{
00298     result.reserve(size());
00299     for (unsigned int ci = beginClusterI() ; ci!= endClusterI();++ci){
00300       SiStripRegionalClusterRef clusterRef = edm::makeRefToLazyGetter(regionalHandle(),ci);
00301       if (isMasked(*clusterRef)) continue;
00302       if (accept(clusterRef))
00303         buildSimpleRecHit( clusterRef, ts,result);
00304       else LogDebug("TkStripMeasurementDet")<<"skipping this reg str from last iteration on"<<rawId()<<" key: "<<clusterRef.key();
00305     }
00306   }
00307 }
00308 
00309 
00310 
00311 bool
00312 TkStripMeasurementDet::testStrips(float utraj, float uerr) const {
00313     int16_t start = (int16_t) std::max<float>(utraj - 3.f*uerr, 0);
00314     int16_t end   = (int16_t) std::min<float>(utraj + 3.f*uerr, totalStrips());
00315 
00316     if (start >= end) { // which means either end <=0 or start >= totalStrips_
00317         /* LogDebug("TkStripMeasurementDet") << "Testing module " << id_ <<","<<
00318             " U = " << utraj << " +/- " << uerr << 
00319             "; Range [" << (utraj - 3*uerr) << ", " << (utraj + 3*uerr) << "] " << 
00320             ": YOU'RE COMPLETELY OFF THE MODULE."; */
00321         //return false; 
00322         return true;  // Wolfgang thinks this way is better
00323                       // and solves some problems with grouped ckf
00324     } 
00325 
00326     typedef std::vector<BadStripBlock>::const_iterator BSBIT;
00327 
00328     int16_t bad = 0, largestBadBlock = 0;
00329     for (BSBIT bsbc = badStripBlocks().begin(), bsbe = badStripBlocks().end(); bsbc != bsbe; ++bsbc) {
00330         if (bsbc->last  < start) continue;
00331         if (bsbc->first > end)   break;
00332         int16_t thisBad = std::min(bsbc->last, end) - std::max(bsbc->first, start);
00333         if (thisBad > largestBadBlock) largestBadBlock = thisBad;
00334         bad += thisBad;
00335     }
00336 
00337     bool ok = (bad < (end-start)) && 
00338       (uint16_t(bad) <= badStripCuts().maxBad) && 
00339       (uint16_t(largestBadBlock) <= badStripCuts().maxConsecutiveBad);
00340 
00341 //    if (bad) {   
00342 //       edm::LogWarning("TkStripMeasurementDet") << "Testing module " << id_ <<" (subdet: "<< SiStripDetId(id_).subdetId() << ")" <<
00343 //            " U = " << utraj << " +/- " << uerr << 
00344 //            "; Range [" << (utraj - 3*uerr) << ", " << (utraj + 3*uerr) << "] " << 
00345 //            "= [" << start << "," << end << "]" <<
00346 //            " total strips:" << (end-start) << ", good:" << (end-start-bad) << ", bad:" << bad << ", largestBadBlock: " << largestBadBlock << 
00347 //            ". " << (ok ? "OK" : "NO"); 
00348 //    }
00349     return ok;
00350 }
00351