CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/RecoTracker/MeasurementDet/src/MeasurementTracker.cc

Go to the documentation of this file.
00001 #include "RecoTracker/MeasurementDet/interface/MeasurementTracker.h"
00002 
00003 #include "FWCore/Framework/interface/ESHandle.h"
00004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00005 
00006 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00007 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
00008 #include "Geometry/CommonDetUnit/interface/GeomDetUnit.h"
00009 #include "Geometry/TrackerGeometryBuilder/interface/GluedGeomDet.h"
00010 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00011 
00012 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00013 
00014 #include "DataFormats/DetId/interface/DetIdCollection.h"
00015 
00016 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
00017 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
00018 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
00019 #include "DataFormats/SiStripCluster/interface/SiStripClusterCollection.h"
00020 
00021 #include "TrackingTools/MeasurementDet/interface/MeasurementDetException.h"
00022 
00023 #include "RecoLocalTracker/ClusterParameterEstimator/interface/PixelClusterParameterEstimator.h"
00024 #include "RecoLocalTracker/Records/interface/TrackerCPERecord.h"
00025 #include "RecoLocalTracker/SiStripRecHitConverter/interface/SiStripRecHitMatcher.h"
00026 #include "RecoLocalTracker/SiStripRecHitConverter/interface/StripCPE.h"  
00027 
00028 #include "RecoTracker/TkDetLayers/interface/GeometricSearchTracker.h"
00029 #include "RecoTracker/Record/interface/TrackerRecoGeometryRecord.h"
00030 #include "RecoTracker/MeasurementDet/interface/TkStripMeasurementDet.h"
00031 #include "RecoTracker/MeasurementDet/interface/TkPixelMeasurementDet.h"
00032 #include "RecoTracker/MeasurementDet/interface/TkGluedMeasurementDet.h"
00033 
00034 #include "CondFormats/SiStripObjects/interface/SiStripNoises.h"
00035 #include "CondFormats/DataRecord/interface/SiStripNoisesRcd.h"
00036 
00037 #include "FWCore/ServiceRegistry/interface/Service.h"
00038 #include "FWCore/Services/interface/UpdaterService.h"
00039 
00040 #include <iostream>
00041 #include <typeinfo>
00042 #include <map>
00043 #include <algorithm>
00044 
00045 
00046 // here just while testing
00047 
00048 #if defined(__GNUC__) && (__GNUC__ == 4) && (__GNUC_MINOR__ > 4)
00049 #include <x86intrin.h>
00050 
00051 #else
00052 
00053 #ifdef __SSE2__
00054 #include <mmintrin.h>
00055 #include <emmintrin.h>
00056 #endif
00057 #ifdef __SSE3__
00058 #include <pmmintrin.h>
00059 #endif
00060 #ifdef __SSE4_1__
00061 #include <smmintrin.h>
00062 #endif
00063 
00064 #endif
00065 
00066 
00067 //
00068 
00069 using namespace std;
00070 
00071 namespace {
00072 
00073   struct CmpTKD {
00074     bool operator()(MeasurementDet* rh, MeasurementDet* lh) {
00075       return rh->geomDet().geographicalId().rawId() < lh->geomDet().geographicalId().rawId();
00076     }
00077   };
00078 
00079   template<typename TKD>
00080   void sortTKD( std::vector<TKD*> & det) {
00081     std::sort(det.begin(),det.end(),CmpTKD());
00082   }
00083 }
00084 
00085 
00086 MeasurementTracker::MeasurementTracker(const edm::ParameterSet&              conf,
00087                                        const PixelClusterParameterEstimator* pixelCPE,
00088                                        const StripClusterParameterEstimator* stripCPE,
00089                                        const SiStripRecHitMatcher*  hitMatcher,
00090                                        const TrackerGeometry*  trackerGeom,
00091                                        const GeometricSearchTracker* geometricSearchTracker,
00092                                        const SiStripQuality *stripQuality,
00093                                        int   stripQualityFlags,
00094                                        int   stripQualityDebugFlags,
00095                                        const SiPixelQuality *pixelQuality,
00096                                        const SiPixelFedCabling *pixelCabling,
00097                                        int   pixelQualityFlags,
00098                                        int   pixelQualityDebugFlags,
00099                                        bool isRegional) :
00100   pset_(conf),
00101   name_(conf.getParameter<std::string>("ComponentName")),
00102   thePixelCPE(pixelCPE),theStripCPE(stripCPE),theHitMatcher(hitMatcher),
00103   theTrackerGeom(trackerGeom),theGeometricSearchTracker(geometricSearchTracker),
00104   theInactivePixelDetectorLabels(conf.getParameter<std::vector<edm::InputTag> >("inactivePixelDetectorLabels")),
00105   theInactiveStripDetectorLabels(conf.getParameter<std::vector<edm::InputTag> >("inactiveStripDetectorLabels")),
00106   isRegional_(isRegional)
00107 {
00108   this->initialize();
00109   this->initializeStripStatus(stripQuality, stripQualityFlags, stripQualityDebugFlags);
00110   this->initializePixelStatus(pixelQuality, pixelCabling, pixelQualityFlags, pixelQualityDebugFlags);
00111   //the measurement tracking is set to skip clusters, the other option is set from outside
00112   selfUpdateSkipClusters_=conf.exists("skipClusters");
00113   if (selfUpdateSkipClusters_)
00114     {
00115       edm::InputTag skip=conf.getParameter<edm::InputTag>("skipClusters");
00116       if (skip==edm::InputTag("")) selfUpdateSkipClusters_=false;
00117     }
00118   
00119   LogDebug("MeasurementTracker")<<"skipping clusters: "<<selfUpdateSkipClusters_;
00120 }
00121 
00122 MeasurementTracker::~MeasurementTracker()
00123 {
00124   for(vector<TkPixelMeasurementDet*>::const_iterator it=thePixelDets.begin(); it!=thePixelDets.end(); ++it){
00125     delete *it;
00126   }
00127 
00128   for(vector<TkStripMeasurementDet*>::const_iterator it=theStripDets.begin(); it!=theStripDets.end(); ++it){
00129     delete *it;
00130   }
00131 
00132   for(vector<TkGluedMeasurementDet*>::const_iterator it=theGluedDets.begin(); it!=theGluedDets.end(); ++it){
00133     delete *it;
00134   }
00135 }
00136 
00137 
00138 void MeasurementTracker::initialize() const
00139 {  
00140   addPixelDets( theTrackerGeom->detsPXB());
00141   addPixelDets( theTrackerGeom->detsPXF());
00142   addStripDets( theTrackerGeom->detsTIB());
00143   addStripDets( theTrackerGeom->detsTID());
00144   addStripDets( theTrackerGeom->detsTOB());
00145   addStripDets( theTrackerGeom->detsTEC());  
00146 
00147   sortTKD(theStripDets);
00148   sortTKD(thePixelDets);
00149   sortTKD(theGluedDets);
00150 
00151 }
00152 
00153 
00154 void MeasurementTracker::addPixelDets( const TrackingGeometry::DetContainer& dets) const
00155 {
00156   for (TrackerGeometry::DetContainer::const_iterator gd=dets.begin();
00157        gd != dets.end(); gd++) {
00158     addPixelDet(*gd, thePixelCPE);
00159   }  
00160 }
00161 
00162 void MeasurementTracker::addStripDets( const TrackingGeometry::DetContainer& dets) const
00163 {
00164   for (TrackerGeometry::DetContainer::const_iterator gd=dets.begin();
00165        gd != dets.end(); gd++) {
00166 
00167     const GeomDetUnit* gdu = dynamic_cast<const GeomDetUnit*>(*gd);
00168 
00169     //    StripSubdetector stripId( (**gd).geographicalId());
00170     //     bool isDetUnit( gdu != 0);
00171     //     cout << "StripSubdetector glued? " << stripId.glued() 
00172     //   << " is DetUnit? " << isDetUnit << endl;
00173 
00174     if (gdu != 0) {
00175       addStripDet(*gd, theStripCPE);
00176     }
00177     else {
00178       const GluedGeomDet* gluedDet = dynamic_cast<const GluedGeomDet*>(*gd);
00179       if (gluedDet == 0) {
00180         throw MeasurementDetException("MeasurementTracker ERROR: GeomDet neither DetUnit nor GluedDet");
00181       }
00182       addGluedDet(gluedDet, theHitMatcher);
00183     }  
00184   }
00185 }
00186 
00187 void MeasurementTracker::addStripDet( const GeomDet* gd,
00188                                       const StripClusterParameterEstimator* cpe) const
00189 {
00190   try {
00191     TkStripMeasurementDet* det = new TkStripMeasurementDet( gd, cpe,isRegional_);
00192     theStripDets.push_back(det);
00193     theDetMap[gd->geographicalId()] = det;
00194   }
00195   catch(MeasurementDetException& err){
00196     edm::LogError("MeasurementDet") << "Oops, got a MeasurementDetException: " << err.what() ;
00197   }
00198 }
00199 
00200 void MeasurementTracker::addPixelDet( const GeomDet* gd,
00201                                       const PixelClusterParameterEstimator* cpe) const
00202 {
00203   TkPixelMeasurementDet* det = new TkPixelMeasurementDet( gd, cpe);
00204   thePixelDets.push_back(det);
00205   theDetMap[gd->geographicalId()] = det;
00206 }
00207 
00208 void MeasurementTracker::addGluedDet( const GluedGeomDet* gd,
00209                                       const SiStripRecHitMatcher* matcher) const
00210 {
00211   const MeasurementDet* monoDet = idToDet( gd->monoDet()->geographicalId());
00212   if (monoDet == 0) {
00213     addStripDet(gd->monoDet(), theStripCPE);  // in case glued det comes before components
00214     monoDet = idToDet( gd->monoDet()->geographicalId());
00215   }
00216 
00217   const MeasurementDet* stereoDet = idToDet( gd->stereoDet()->geographicalId());
00218   if (stereoDet == 0) {
00219     addStripDet(gd->stereoDet(), theStripCPE);  // in case glued det comes before components
00220     stereoDet = idToDet( gd->stereoDet()->geographicalId());
00221   }
00222 
00223   if (monoDet == 0 || stereoDet == 0) {
00224     edm::LogError("MeasurementDet") << "MeasurementTracker ERROR: GluedDet components not found as MeasurementDets ";
00225     throw MeasurementDetException("MeasurementTracker ERROR: GluedDet components not found as MeasurementDets");
00226   }
00227 
00228   TkGluedMeasurementDet* det = new TkGluedMeasurementDet( gd, theHitMatcher,
00229                                                           monoDet, stereoDet);
00230   theGluedDets.push_back( det);
00231   theDetMap[gd->geographicalId()] = det;
00232 }
00233 
00234 void MeasurementTracker::update( const edm::Event& event) const
00235 {
00236   updatePixels(event);
00237   updateStrips(event);
00238   
00239   /*
00240   for (std::vector<TkStripMeasurementDet*>::const_iterator i=theStripDets.begin();
00241        i!=theStripDets.end(); i++) {
00242     if( (*i)->isEmpty()){
00243       std::cout << "stripDet id, #hits: " 
00244                 <<  (*i)->geomDet().geographicalId().rawId() << " , "
00245                 << 0 << std::endl;
00246     }else{
00247       std::cout << "stripDet id, #hits: " 
00248                 <<  (*i)->geomDet().geographicalId().rawId() << " , "
00249                 << (*i)->size() << std::endl;
00250     }
00251   }
00252   */
00253 }
00254 void MeasurementTracker::setClusterToSkip(const edm::InputTag & cluster, const edm::Event& event) const
00255 {
00256   LogDebug("MeasurementTracker")<<"setClusterToSkip";
00257   if (selfUpdateSkipClusters_)
00258     edm::LogError("MeasurementTracker")<<"this mode of operation is not supported, either the measurement tracker is set to skip clusters, or is being told to skip clusters. not both";
00259 
00260   edm::Handle<edmNew::DetSetVector<SiPixelClusterRefNew> > pixelClusterRefs;
00261   event.getByLabel(cluster,pixelClusterRefs);
00262   for (std::vector<TkPixelMeasurementDet*>::const_iterator i=thePixelDets.begin();
00263        i!=thePixelDets.end(); i++) {
00264     edmNew::DetSetVector<SiPixelClusterRefNew>::const_iterator f=pixelClusterRefs->find((**i).geomDet().geographicalId().rawId());
00265     if (f!=pixelClusterRefs->end())
00266       (**i).setClusterToSkip(f->begin(),f->end());
00267     else
00268       (**i).unset();
00269   }
00270 
00271   edm::Handle<edmNew::DetSetVector<TkStripMeasurementDet::SiStripClusterRef> > stripClusterRefs;
00272   event.getByLabel(cluster,stripClusterRefs);
00273   for (std::vector<TkStripMeasurementDet*>::const_iterator i=theStripDets.begin();
00274        i!=theStripDets.end(); i++){
00275     edmNew::DetSetVector<TkStripMeasurementDet::SiStripClusterRef>::const_iterator f=stripClusterRefs->find((**i).geomDet().geographicalId().rawId());
00276     if (f!=stripClusterRefs->end())
00277       (**i).setClusterToSkip(f->begin(),f->end());
00278     else
00279       (**i).unset();
00280   }
00281 }
00282 
00283 void MeasurementTracker::unsetClusterToSkip() const {
00284   LogDebug("MeasurementTracker")<<"unsetClusterToSkip";
00285   if (selfUpdateSkipClusters_)
00286     edm::LogError("MeasurementTracker")<<"this mode of operation is not supported, either the measurement tracker is set to skip clusters, or is being told to skip clusters. not both";
00287 
00288   for (std::vector<TkPixelMeasurementDet*>::const_iterator i=thePixelDets.begin();
00289        i!=thePixelDets.end(); i++) {
00290     (**i).unset();
00291   }
00292   for (std::vector<TkStripMeasurementDet*>::const_iterator i=theStripDets.begin();
00293        i!=theStripDets.end(); i++){
00294     (**i).unset();
00295   }
00296 }
00297 
00298 void MeasurementTracker::updatePixels( const edm::Event& event) const
00299 {
00300   // avoid to update twice from the same event
00301   if (!edm::Service<UpdaterService>()->checkOnce("MeasurementTracker::updatePixels::"+name_)) return;
00302 
00303   typedef edmNew::DetSet<SiPixelCluster> PixelDetSet;
00304 
00305   bool switchOffPixelsIfEmpty = (!pset_.existsAs<bool>("switchOffPixelsIfEmpty")) ||
00306                                 (pset_.getParameter<bool>("switchOffPixelsIfEmpty"));
00307 
00308   std::vector<uint32_t> rawInactiveDetIds; 
00309   if (!theInactivePixelDetectorLabels.empty()) {
00310     edm::Handle<DetIdCollection> detIds;
00311     for (std::vector<edm::InputTag>::const_iterator itt = theInactivePixelDetectorLabels.begin(), edt = theInactivePixelDetectorLabels.end(); 
00312             itt != edt; ++itt) {
00313         event.getByLabel(*itt, detIds);
00314         rawInactiveDetIds.insert(rawInactiveDetIds.end(), detIds->begin(), detIds->end());
00315     }
00316     if (!rawInactiveDetIds.empty()) std::sort(rawInactiveDetIds.begin(), rawInactiveDetIds.end());
00317   }
00318   // Pixel Clusters
00319   std::string pixelClusterProducer = pset_.getParameter<std::string>("pixelClusterProducer");
00320   if( pixelClusterProducer.empty() ) { //clusters have not been produced
00321     for (std::vector<TkPixelMeasurementDet*>::const_iterator i=thePixelDets.begin();
00322          i!=thePixelDets.end(); i++) {
00323       if (switchOffPixelsIfEmpty) {
00324         (**i).setActiveThisEvent(false);
00325       }else{
00326         (**i).setEmpty();
00327       }
00328     }
00329   }else{  
00330 
00331     edm::Handle<edmNew::DetSetVector<SiPixelClusterRefNew> > pixelClusterRefs;
00332     if (selfUpdateSkipClusters_){
00333       //and get the collection of pixel ref to skip
00334       event.getByLabel(pset_.getParameter<edm::InputTag>("skipClusters"),pixelClusterRefs);
00335       LogDebug("MeasurementTracker")<<"getting pxl refs to skip";
00336       if (pixelClusterRefs.failedToGet())edm::LogError("MeasurementTracker")<<"not getting the pixel clusters to skip";
00337     }
00338 
00339     edm::Handle<edmNew::DetSetVector<SiPixelCluster> > pixelClusters;
00340     event.getByLabel(pixelClusterProducer, pixelClusters);
00341     const  edmNew::DetSetVector<SiPixelCluster>* pixelCollection = pixelClusters.product();
00342    
00343     if (switchOffPixelsIfEmpty && pixelCollection->empty()) {
00344         for (std::vector<TkPixelMeasurementDet*>::const_iterator i=thePixelDets.begin();
00345              i!=thePixelDets.end(); i++) {
00346               (**i).setActiveThisEvent(false);
00347         }
00348     } else { 
00349 
00350     for (std::vector<TkPixelMeasurementDet*>::const_iterator i=thePixelDets.begin();
00351          i!=thePixelDets.end(); i++) {
00352 
00353       // foreach det get cluster range
00354       unsigned int id = (**i).geomDet().geographicalId().rawId();
00355       if (!rawInactiveDetIds.empty() && std::binary_search(rawInactiveDetIds.begin(), rawInactiveDetIds.end(), id)) {
00356         (**i).setActiveThisEvent(false); continue;
00357       }
00358       //FIXME
00359       //fill the set with what needs to be skipped
00360       edmNew::DetSetVector<SiPixelCluster>::const_iterator it = pixelCollection->find( id );
00361       if ( it != pixelCollection->end() ){            
00362         // push cluster range in det
00363         (**i).update( *it, pixelClusters, id );
00364         if (selfUpdateSkipClusters_){
00365           edmNew::DetSetVector<SiPixelClusterRefNew>::const_iterator f=pixelClusterRefs->find(id);
00366           if (f!=pixelClusterRefs->end())
00367             (**i).setClusterToSkip(f->begin(),f->end());
00368         }
00369       }else{
00370         (**i).setEmpty();
00371       }
00372     }
00373     }
00374   }
00375 
00376 
00377   
00378 }
00379 
00380 void MeasurementTracker::getInactiveStrips(const edm::Event& event,
00381                                            std::vector<uint32_t> & rawInactiveDetIds) const {
00382   if (!theInactiveStripDetectorLabels.empty()) {
00383     edm::Handle<DetIdCollection> detIds;
00384     for (std::vector<edm::InputTag>::const_iterator itt = theInactiveStripDetectorLabels.begin(), edt = theInactiveStripDetectorLabels.end();
00385          itt != edt; ++itt) {
00386       event.getByLabel(*itt, detIds);
00387       rawInactiveDetIds.insert(rawInactiveDetIds.end(), detIds->begin(), detIds->end());
00388     }
00389     if (!rawInactiveDetIds.empty()) std::sort(rawInactiveDetIds.begin(), rawInactiveDetIds.end());
00390   }
00391 }
00392 
00393 void MeasurementTracker::updateStrips( const edm::Event& event) const
00394 {
00395   // avoid to update twice from the same event
00396   if (!edm::Service<UpdaterService>()->checkOnce("MeasurementTracker::updateStrips::"+name_)) return;
00397 
00398   typedef edmNew::DetSet<SiStripCluster>   StripDetSet;
00399 
00400   std::vector<uint32_t> rawInactiveDetIds;
00401   getInactiveStrips(event,rawInactiveDetIds);
00402 
00403   // Strip Clusters
00404   std::string stripClusterProducer = pset_.getParameter<std::string>("stripClusterProducer");
00405   //first clear all of them
00406 
00407   {
00408     std::vector<TkStripMeasurementDet*>::const_iterator end = theStripDets.end()-200;
00409     for (std::vector<TkStripMeasurementDet*>::const_iterator i=theStripDets.begin();
00410          i!=end; i++) {
00411       (**i).setEmpty();
00412 #ifdef __SSE2__
00413       _mm_prefetch(((char *)(*(i+200))),_MM_HINT_T0); 
00414 #endif
00415     }
00416    for (std::vector<TkStripMeasurementDet*>::const_iterator i=end;
00417          i!=theStripDets.end(); i++)
00418       (**i).setEmpty();
00419   }
00420   if( !stripClusterProducer.compare("") ) { //clusters have not been produced
00421   }else{
00422     //=========  actually load cluster =============
00423     if(!isRegional_){
00424         edm::Handle<edmNew::DetSetVector<TkStripMeasurementDet::SiStripClusterRef> > stripClusterRefs;
00425         if (selfUpdateSkipClusters_){
00426           //and get the collection of pixel ref to skip
00427           LogDebug("MeasurementTracker")<<"getting strp refs to skip";
00428           event.getByLabel(pset_.getParameter<edm::InputTag>("skipClusters"),stripClusterRefs);
00429           if (stripClusterRefs.failedToGet())edm::LogError("MeasurementTracker")<<"not getting the strip clusters to skip";
00430         }
00431 
00432       edm::Handle<edmNew::DetSetVector<SiStripCluster> > clusterHandle;
00433       event.getByLabel(stripClusterProducer, clusterHandle);
00434       const edmNew::DetSetVector<SiStripCluster>* clusterCollection = clusterHandle.product();
00435 
00436       std::vector<TkStripMeasurementDet*>::const_iterator i=theStripDets.begin();
00437       std::vector<TkStripMeasurementDet*>::const_iterator endDet=theStripDets.end();
00438       edmNew::DetSetVector<SiStripCluster>::const_iterator it = (*clusterCollection).begin();
00439       edmNew::DetSetVector<SiStripCluster>::const_iterator endColl = (*clusterCollection).end();
00440       // cluster and det and in order (both) and unique so let's use set intersection
00441       for (;it!=endColl; ++it) {
00442         StripDetSet detSet = *it;
00443         unsigned int id = detSet.id();
00444         while ( id != (**i).rawId()) { // eventually change to lower_range
00445           ++i;
00446           if (i==endDet) throw "we have a problem!!!!";
00447         }
00448 
00449         if (!rawInactiveDetIds.empty() && std::binary_search(rawInactiveDetIds.begin(), rawInactiveDetIds.end(), id)) {
00450             (**i).setActiveThisEvent(false); continue;
00451         }
00452         // push cluster range in det
00453 
00454         (**i).update( detSet, clusterHandle, id );
00455         if (selfUpdateSkipClusters_){
00456           edmNew::DetSetVector<TkStripMeasurementDet::SiStripClusterRef>::const_iterator f=stripClusterRefs->find(id);
00457           if (f!=stripClusterRefs->end())
00458             (**i).setClusterToSkip(f->begin(),f->end());
00459         }
00460       }
00461     }else{
00462 
00463       edm::Handle<edmNew::DetSetVector<TkStripMeasurementDet::SiStripRegionalClusterRef> > stripClusterRefs;
00464       if(selfUpdateSkipClusters_){
00465         LogDebug("MeasurementTracker")<<"getting reg strp refs to skip";
00466         event.getByLabel(pset_.getParameter<edm::InputTag>("skipClusters"),stripClusterRefs);
00467         if (stripClusterRefs.failedToGet())edm::LogError("MeasurementTracker")<<"not getting the strip clusters to skip";
00468       }
00469 
00470       //then set the not-empty ones only
00471       edm::Handle<edm::RefGetter<SiStripCluster> > refClusterHandle;
00472       event.getByLabel(stripClusterProducer, refClusterHandle);
00473       
00474       std::string lazyGetter = pset_.getParameter<std::string>("stripLazyGetterProducer");
00475       edm::Handle<edm::LazyGetter<SiStripCluster> > lazyClusterHandle;
00476       event.getByLabel(lazyGetter,lazyClusterHandle);
00477 
00478       uint32_t tmpId=0;
00479       vector<SiStripCluster>::const_iterator beginIterator;
00480       edm::RefGetter<SiStripCluster>::const_iterator iregion = refClusterHandle->begin();
00481       for(;iregion!=refClusterHandle->end();++iregion) {
00482         const edm::RegionIndex<SiStripCluster>& region = *iregion;
00483         vector<SiStripCluster>::const_iterator icluster = region.begin();
00484         const vector<SiStripCluster>::const_iterator endIterator = region.end();
00485         tmpId = icluster->geographicalId();
00486         beginIterator = icluster;
00487 
00488         //std::cout << "== tmpId ad inizio loop dentro region: " << tmpId << std::endl;
00489 
00490         for (;icluster!=endIterator;icluster++) {
00491           //std::cout << "===== cluster id,pos " 
00492           //  << icluster->geographicalId() << " , " << icluster->barycenter()
00493           //  << std::endl;
00494           //std::cout << "=====making ref in recHits() " << std::endl;
00495           if( icluster->geographicalId() != tmpId){ 
00496             //std::cout << "geo!=tmpId" << std::endl;
00497             //we should find a way to avoid this casting. it is slow
00498             //create also another map for TkStripMeasurementDet ??
00499 
00500             // the following castings are really ugly. To be corrected ASAP
00501             const TkStripMeasurementDet* theConcreteDet = 
00502               dynamic_cast<const TkStripMeasurementDet*>(idToDet(DetId(tmpId)));
00503             
00504             if(theConcreteDet == 0)
00505               throw MeasurementDetException("failed casting to TkStripMeasurementDet*");            
00506             TkStripMeasurementDet*  theConcreteDetUpdatable = 
00507               const_cast<TkStripMeasurementDet*>(theConcreteDet);
00508             theConcreteDetUpdatable->update(beginIterator,icluster,lazyClusterHandle,tmpId);
00509             if (selfUpdateSkipClusters_){
00510               edmNew::DetSetVector<TkStripMeasurementDet::SiStripRegionalClusterRef>::const_iterator f=stripClusterRefs->find(tmpId);
00511               if (f!=stripClusterRefs->end())
00512                 theConcreteDetUpdatable->setRegionalClustersToSkip(f->begin(),f->end());
00513             }
00514             //cannot we avoid to update the det with detId of itself??
00515 
00516             tmpId = icluster->geographicalId();
00517             beginIterator = icluster;
00518             if( icluster == (endIterator-1)){
00519               const TkStripMeasurementDet* theConcreteDet = 
00520               dynamic_cast<const TkStripMeasurementDet*>(idToDet(DetId(tmpId)));
00521               
00522               if(theConcreteDet == 0)
00523               throw MeasurementDetException("failed casting to TkStripMeasurementDet*");            
00524               TkStripMeasurementDet*  theConcreteDetUpdatable = 
00525               const_cast<TkStripMeasurementDet*>(theConcreteDet);
00526               theConcreteDetUpdatable->update(icluster,endIterator,lazyClusterHandle,tmpId);
00527               if (selfUpdateSkipClusters_){
00528                 edmNew::DetSetVector<TkStripMeasurementDet::SiStripRegionalClusterRef>::const_iterator f=stripClusterRefs->find(tmpId);
00529                 if (f!=stripClusterRefs->end())
00530                   theConcreteDetUpdatable->setRegionalClustersToSkip(f->begin(),f->end());
00531               }
00532             }   
00533           }else if( icluster == (endIterator-1)){          
00534             const TkStripMeasurementDet* theConcreteDet = 
00535               dynamic_cast<const TkStripMeasurementDet*>(idToDet(DetId(tmpId)));
00536             
00537             if(theConcreteDet == 0)
00538               throw MeasurementDetException("failed casting to TkStripMeasurementDet*");            
00539             TkStripMeasurementDet*  theConcreteDetUpdatable = 
00540               const_cast<TkStripMeasurementDet*>(theConcreteDet);
00541             //std::cout << "=== option3. fill det with id,#clust: " << tmpId  << " , " 
00542             //      << iregion->end() - beginIterator << std::endl;
00543             theConcreteDetUpdatable->update(beginIterator,endIterator,lazyClusterHandle,tmpId);  
00544             if (selfUpdateSkipClusters_){
00545               edmNew::DetSetVector<TkStripMeasurementDet::SiStripRegionalClusterRef>::const_iterator f=stripClusterRefs->find(tmpId);
00546               if (f!=stripClusterRefs->end())
00547                 theConcreteDetUpdatable->setRegionalClustersToSkip(f->begin(),f->end());
00548             }    
00549           }
00550         }//end loop cluster in one ragion
00551       }
00552     }//end of block for updating with regional clusters 
00553   }
00554 }
00555 
00556 
00557 
00558 const MeasurementDet* 
00559 MeasurementTracker::idToDet(const DetId& id) const
00560 {
00561   std::map<DetId,MeasurementDet*>::const_iterator it = theDetMap.find(id);
00562   if(it !=theDetMap.end()) {
00563     return it->second;
00564   }else{
00565     //throw exception;
00566   }
00567   
00568   return 0; //to avoid compile warning
00569 }
00570 
00571 void MeasurementTracker::initializeStripStatus(const SiStripQuality *quality, int qualityFlags, int qualityDebugFlags) const {
00572   TkStripMeasurementDet::BadStripCuts badStripCuts[4];
00573   if (qualityFlags & BadStrips) {
00574      edm::ParameterSet cutPset = pset_.getParameter<edm::ParameterSet>("badStripCuts");
00575      badStripCuts[SiStripDetId::TIB-3] = TkStripMeasurementDet::BadStripCuts(cutPset.getParameter<edm::ParameterSet>("TIB"));
00576      badStripCuts[SiStripDetId::TOB-3] = TkStripMeasurementDet::BadStripCuts(cutPset.getParameter<edm::ParameterSet>("TOB"));
00577      badStripCuts[SiStripDetId::TID-3] = TkStripMeasurementDet::BadStripCuts(cutPset.getParameter<edm::ParameterSet>("TID"));
00578      badStripCuts[SiStripDetId::TEC-3] = TkStripMeasurementDet::BadStripCuts(cutPset.getParameter<edm::ParameterSet>("TEC"));
00579   }
00580 
00581   if ((quality != 0) && (qualityFlags != 0))  {
00582     edm::LogInfo("MeasurementTracker") << "qualityFlags = " << qualityFlags;
00583     unsigned int on = 0, tot = 0; 
00584     unsigned int foff = 0, ftot = 0, aoff = 0, atot = 0; 
00585     for (std::vector<TkStripMeasurementDet*>::const_iterator i=theStripDets.begin();
00586          i!=theStripDets.end(); i++) {
00587       uint32_t detid = ((**i).geomDet().geographicalId()).rawId();
00588       if (qualityFlags & BadModules) {
00589           bool isOn = quality->IsModuleUsable(detid);
00590           (*i)->setActive(isOn);
00591           tot++; on += (unsigned int) isOn;
00592           if (qualityDebugFlags & BadModules) {
00593             edm::LogInfo("MeasurementTracker")<< "MeasurementTracker::initializeStripStatus : detid " << detid << " is " << (isOn ?  "on" : "off");
00594           }
00595        } else {
00596           (*i)->setActive(true);
00597        }
00598        // first turn all APVs and fibers ON
00599        (*i)->set128StripStatus(true); 
00600        if (qualityFlags & BadAPVFibers) {
00601           short badApvs   = quality->getBadApvs(detid);
00602           short badFibers = quality->getBadFibers(detid);
00603           for (int j = 0; j < 6; j++) {
00604              atot++;
00605              if (badApvs & (1 << j)) {
00606                 (*i)->set128StripStatus(false, j);
00607                 aoff++;
00608              }
00609           }
00610           for (int j = 0; j < 3; j++) {
00611              ftot++;
00612              if (badFibers & (1 << j)) {
00613                 (*i)->set128StripStatus(false, 2*j);
00614                 (*i)->set128StripStatus(false, 2*j+1);
00615                 foff++;
00616              }
00617           }
00618           (*i)->setMaskBad128StripBlocks((qualityFlags & MaskBad128StripBlocks) != 0);
00619        } 
00620        std::vector<TkStripMeasurementDet::BadStripBlock> &badStrips = (*i)->getBadStripBlocks();
00621        badStrips.clear();
00622        if (qualityFlags & BadStrips) {
00623             SiStripBadStrip::Range range = quality->getRange(detid);
00624             for (SiStripBadStrip::ContainerIterator bit = range.first; bit != range.second; ++bit) {
00625                 badStrips.push_back(quality->decode(*bit));
00626             }
00627             (*i)->setBadStripCuts(badStripCuts[SiStripDetId(detid).subdetId()-3]);
00628        }
00629     }
00630     if (qualityDebugFlags & BadModules) {
00631         edm::LogInfo("MeasurementTracker StripModuleStatus") << 
00632             " Total modules: " << tot << ", active " << on <<", inactive " << (tot - on);
00633     }
00634     if (qualityDebugFlags & BadAPVFibers) {
00635         edm::LogInfo("MeasurementTracker StripAPVStatus") << 
00636             " Total APVs: " << atot << ", active " << (atot-aoff) <<", inactive " << (aoff);
00637         edm::LogInfo("MeasurementTracker StripFiberStatus") << 
00638             " Total Fibers: " << ftot << ", active " << (ftot-foff) <<", inactive " << (foff);
00639     }
00640   } else {
00641     for (std::vector<TkStripMeasurementDet*>::const_iterator i=theStripDets.begin();
00642          i!=theStripDets.end(); i++) {
00643       (*i)->setActive(true);          // module ON
00644       (*i)->set128StripStatus(true);  // all APVs and fibers ON
00645     }
00646   }
00647 }
00648 
00649 void MeasurementTracker::initializePixelStatus(const SiPixelQuality *quality, const SiPixelFedCabling *pixelCabling, int qualityFlags, int qualityDebugFlags) const {
00650   if ((quality != 0) && (qualityFlags != 0))  {
00651     edm::LogInfo("MeasurementTracker") << "qualityFlags = " << qualityFlags;
00652     unsigned int on = 0, tot = 0, badrocs = 0; 
00653     for (std::vector<TkPixelMeasurementDet*>::const_iterator i=thePixelDets.begin();
00654          i!=thePixelDets.end(); i++) {
00655       uint32_t detid = ((**i).geomDet().geographicalId()).rawId();
00656       if (qualityFlags & BadModules) {
00657           bool isOn = quality->IsModuleUsable(detid);
00658           (*i)->setActive(isOn);
00659           tot++; on += (unsigned int) isOn;
00660           if (qualityDebugFlags & BadModules) {
00661             edm::LogInfo("MeasurementTracker")<< "MeasurementTracker::initializePixelStatus : detid " << detid << " is " << (isOn ?  "on" : "off");
00662           }
00663        } else {
00664           (*i)->setActive(true);
00665        }
00666        if ((qualityFlags & BadROCs) && (quality->getBadRocs(detid) != 0)) {
00667           std::vector<LocalPoint> badROCs = quality->getBadRocPositions(detid, *theTrackerGeom, pixelCabling);
00668           badrocs += badROCs.size();
00669           (*i)->setBadRocPositions(badROCs);
00670        } else {
00671           (*i)->clearBadRocPositions();  
00672        }
00673     }
00674     if (qualityDebugFlags & BadModules) {
00675         edm::LogInfo("MeasurementTracker PixelModuleStatus") << 
00676             " Total modules: " << tot << ", active " << on <<", inactive " << (tot - on);
00677     }
00678     if (qualityDebugFlags & BadROCs) {
00679         edm::LogInfo("MeasurementTracker PixelROCStatus") << " Total of bad ROCs: " << badrocs ;
00680     }
00681   } else {
00682     for (std::vector<TkPixelMeasurementDet*>::const_iterator i=thePixelDets.begin();
00683          i!=thePixelDets.end(); i++) {
00684       (*i)->setActive(true);          // module ON
00685     }
00686   }
00687 }
00688