CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/RecoTracker/MeasurementDet/plugins/MeasurementTrackerImpl.cc

Go to the documentation of this file.
00001 #include "MeasurementTrackerImpl.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 #include "DataFormats/Common/interface/ContainerMask.h"
00021 
00022 #include "TrackingTools/MeasurementDet/interface/MeasurementDetException.h"
00023 
00024 #include "RecoLocalTracker/ClusterParameterEstimator/interface/PixelClusterParameterEstimator.h"
00025 #include "RecoLocalTracker/Records/interface/TrackerCPERecord.h"
00026 #include "RecoLocalTracker/SiStripRecHitConverter/interface/SiStripRecHitMatcher.h"
00027 #include "RecoLocalTracker/SiStripRecHitConverter/interface/StripCPE.h"  
00028 
00029 #include "RecoTracker/TkDetLayers/interface/GeometricSearchTracker.h"
00030 #include "RecoTracker/Record/interface/TrackerRecoGeometryRecord.h"
00031 #include "TkStripMeasurementDet.h"
00032 #include "TkPixelMeasurementDet.h"
00033 #include "TkGluedMeasurementDet.h"
00034 
00035 #include "CondFormats/SiStripObjects/interface/SiStripNoises.h"
00036 #include "CondFormats/DataRecord/interface/SiStripNoisesRcd.h"
00037 
00038 #include "FWCore/ServiceRegistry/interface/Service.h"
00039 #include "FWCore/Services/interface/UpdaterService.h"
00040 
00041 #include <iostream>
00042 #include <typeinfo>
00043 #include <map>
00044 #include <algorithm>
00045 
00046 
00047 //
00048 
00049 using namespace std;
00050 
00051 namespace {
00052 
00053   struct CmpTKD {
00054     bool operator()(MeasurementDet const* rh, MeasurementDet const * lh) {
00055       return rh->fastGeomDet().geographicalId().rawId() < lh->fastGeomDet().geographicalId().rawId();
00056     }
00057     bool operator()(MeasurementDet const & rh, MeasurementDet const &  lh) {
00058       return rh.fastGeomDet().geographicalId().rawId() < lh.fastGeomDet().geographicalId().rawId();
00059     }
00060   };
00061 
00062   template<typename TKD>
00063   void sortTKD( std::vector<TKD*> & det) {
00064     std::sort(det.begin(),det.end(),CmpTKD());
00065   }
00066   template<typename TKD>
00067   void sortTKD( std::vector<TKD> & det) {
00068     std::sort(det.begin(),det.end(),CmpTKD());
00069   }
00070 
00071 }
00072 
00073 
00074 MeasurementTrackerImpl::MeasurementTrackerImpl(const edm::ParameterSet&              conf,
00075                                        const PixelClusterParameterEstimator* pixelCPE,
00076                                        const StripClusterParameterEstimator* stripCPE,
00077                                        const SiStripRecHitMatcher*  hitMatcher,
00078                                        const TrackerGeometry*  trackerGeom,
00079                                        const GeometricSearchTracker* geometricSearchTracker,
00080                                        const SiStripQuality *stripQuality,
00081                                        int   stripQualityFlags,
00082                                        int   stripQualityDebugFlags,
00083                                        const SiPixelQuality *pixelQuality,
00084                                        const SiPixelFedCabling *pixelCabling,
00085                                        int   pixelQualityFlags,
00086                                        int   pixelQualityDebugFlags,
00087                                        bool isRegional) :
00088   MeasurementTracker(trackerGeom,geometricSearchTracker),
00089   pset_(conf),
00090   name_(conf.getParameter<std::string>("ComponentName")),
00091   theStDets(hitMatcher,stripCPE,isRegional),
00092   thePixelCPE(pixelCPE),
00093   theInactivePixelDetectorLabels(conf.getParameter<std::vector<edm::InputTag> >("inactivePixelDetectorLabels")),
00094   theInactiveStripDetectorLabels(conf.getParameter<std::vector<edm::InputTag> >("inactiveStripDetectorLabels"))
00095 {
00096   this->initialize();
00097   this->initializeStripStatus(stripQuality, stripQualityFlags, stripQualityDebugFlags);
00098   this->initializePixelStatus(pixelQuality, pixelCabling, pixelQualityFlags, pixelQualityDebugFlags);
00099   //the measurement tracking is set to skip clusters, the other option is set from outside
00100   selfUpdateSkipClusters_=conf.exists("skipClusters");
00101   if (selfUpdateSkipClusters_)
00102     {
00103       edm::InputTag skip=conf.getParameter<edm::InputTag>("skipClusters");
00104       if (skip==edm::InputTag("")) selfUpdateSkipClusters_=false;
00105     }
00106 
00107 
00108   LogDebug("MeasurementTracker")<<"skipping clusters: "<<selfUpdateSkipClusters_;
00109 }
00110 
00111 MeasurementTrackerImpl::~MeasurementTrackerImpl()
00112 {
00113   for(vector<TkPixelMeasurementDet*>::const_iterator it=thePixelDets.begin(); it!=thePixelDets.end(); ++it){
00114     delete *it;
00115   }
00116   
00117 }
00118 
00119 
00120 void MeasurementTrackerImpl::initialize()
00121 {  
00122   addPixelDets( theTrackerGeom->detsPXB());
00123   addPixelDets( theTrackerGeom->detsPXF());
00124 
00125   addStripDets( theTrackerGeom->detsTIB());
00126   addStripDets( theTrackerGeom->detsTID());
00127   addStripDets( theTrackerGeom->detsTOB());
00128   addStripDets( theTrackerGeom->detsTEC());  
00129 
00130   // fist all stripdets
00131   sortTKD(theStripDets);
00132   theStDets.init(theStripDets);
00133   for (unsigned int i=0; i!=theStripDets.size(); ++i)
00134     theDetMap[theStDets.id(i)] = &theStripDets[i];
00135   
00136   // now the glued dets
00137   sortTKD(theGluedDets);
00138   for (unsigned int i=0; i!=theGluedDets.size(); ++i)
00139     initGluedDet(theGluedDets[i]);
00140 
00141   sortTKD(thePixelDets);
00142 
00143 
00144 }
00145 
00146 
00147 void MeasurementTrackerImpl::addPixelDets( const TrackingGeometry::DetContainer& dets)
00148 {
00149   for (TrackerGeometry::DetContainer::const_iterator gd=dets.begin();
00150        gd != dets.end(); gd++) {
00151     addPixelDet(*gd, thePixelCPE);
00152   }  
00153 }
00154 
00155 void MeasurementTrackerImpl::addStripDets( const TrackingGeometry::DetContainer& dets)
00156 {
00157   for (TrackerGeometry::DetContainer::const_iterator gd=dets.begin();
00158        gd != dets.end(); gd++) {
00159 
00160     const GeomDetUnit* gdu = dynamic_cast<const GeomDetUnit*>(*gd);
00161 
00162     //    StripSubdetector stripId( (**gd).geographicalId());
00163     //     bool isDetUnit( gdu != 0);
00164     //     cout << "StripSubdetector glued? " << stripId.glued() 
00165     //   << " is DetUnit? " << isDetUnit << endl;
00166 
00167     if (gdu != 0) {
00168       addStripDet(*gd);
00169     }
00170     else {
00171       const GluedGeomDet* gluedDet = dynamic_cast<const GluedGeomDet*>(*gd);
00172       if (gluedDet == 0) {
00173         throw MeasurementDetException("MeasurementTracker ERROR: GeomDet neither DetUnit nor GluedDet");
00174       }
00175       addGluedDet(gluedDet);
00176     }  
00177   }
00178 }
00179 
00180 void MeasurementTrackerImpl::addStripDet( const GeomDet* gd)
00181 {
00182   try {
00183     theStripDets.push_back(TkStripMeasurementDet( gd, theStDets));
00184   }
00185   catch(MeasurementDetException& err){
00186     edm::LogError("MeasurementDet") << "Oops, got a MeasurementDetException: " << err.what() ;
00187   }
00188 }
00189 
00190 void MeasurementTrackerImpl::addPixelDet( const GeomDet* gd,
00191                                       const PixelClusterParameterEstimator* cpe)
00192 {
00193   TkPixelMeasurementDet* det = new TkPixelMeasurementDet( gd, cpe);
00194   thePixelDets.push_back(det);
00195   det->setClusterToSkip(&thePixelsToSkip);
00196   theDetMap[gd->geographicalId()] = det;
00197 }
00198 
00199 void MeasurementTrackerImpl::addGluedDet( const GluedGeomDet* gd)
00200 {
00201   theGluedDets.push_back(TkGluedMeasurementDet( gd, theStDets.matcher(), theStDets.stripCPE()));
00202 }
00203 
00204 void MeasurementTrackerImpl::initGluedDet( TkGluedMeasurementDet & det)
00205 {
00206   const GluedGeomDet& gd = det.specificGeomDet();
00207   const MeasurementDet* monoDet = findDet( gd.monoDet()->geographicalId());
00208   const MeasurementDet* stereoDet = findDet( gd.stereoDet()->geographicalId());
00209   if (monoDet == 0 || stereoDet == 0) {
00210     edm::LogError("MeasurementDet") << "MeasurementTracker ERROR: GluedDet components not found as MeasurementDets ";
00211     throw MeasurementDetException("MeasurementTracker ERROR: GluedDet components not found as MeasurementDets");
00212   }
00213   det.init(monoDet,stereoDet);
00214   theDetMap[gd.geographicalId()] = &det;
00215 }
00216 
00217 
00218 void MeasurementTrackerImpl::update( const edm::Event& event) const
00219 {
00220   updatePixels(event);
00221   updateStrips(event);
00222   
00223   /*
00224   for (std::vector<TkStripMeasurementDet>::const_iterator i=theStripDets.begin();
00225        i!=theStripDets.end(); i++) {
00226     if( (*i).isEmpty()){
00227       std::cout << "stripDet id, #hits: " 
00228                 <<  (*i).geomDet().geographicalId().rawId() << " , "
00229                 << 0 << std::endl;
00230     }else{
00231       std::cout << "stripDet id, #hits: " 
00232                 <<  (*i).geomDet().geographicalId().rawId() << " , "
00233                 << (*i).size() << " " << (*i).detSet().size() std::endl;
00234     }
00235   }
00236   */
00237 }
00238 void MeasurementTrackerImpl::setClusterToSkip(const edm::InputTag & cluster, const edm::Event& event) const
00239 {
00240   //method called by user of the measurement tracker to tell what needs to be skiped from the event.
00241   //there it is incompatible with a configuration in which the measurement tracker already knows what to skip
00242   // i.e selfUpdateSkipClusters_=True
00243 
00244   LogDebug("MeasurementTracker")<<"setClusterToSkip";
00245   if (selfUpdateSkipClusters_)
00246     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";
00247 
00248   edm::Handle<edm::ContainerMask<edmNew::DetSetVector<SiPixelCluster> > > pixelClusterMask;
00249   event.getByLabel(cluster,pixelClusterMask);
00250 
00251 
00252   thePixelsToSkip.resize(pixelClusterMask->size());
00253   pixelClusterMask->copyMaskTo(thePixelsToSkip);
00254     
00255   edm::Handle<edm::ContainerMask<edmNew::DetSetVector<SiStripCluster> > > stripClusterMask;
00256   event.getByLabel(cluster,stripClusterMask);
00257 
00258   theStDets.theStripsToSkip.resize(stripClusterMask->size());
00259   stripClusterMask->copyMaskTo(theStDets.theStripsToSkip);
00260   
00261 }
00262 
00263 void MeasurementTrackerImpl::unsetClusterToSkip() const {
00264   //method called by user of the measurement tracker to tell what needs to be skiped from the event.
00265   //there it is incompatible with a configuration in which the measurement tracker already knows what to skip
00266   // i.e selfUpdateSkipClusters_=True
00267   
00268   LogDebug("MeasurementTracker")<<"unsetClusterToSkip";
00269   if (selfUpdateSkipClusters_)
00270     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";
00271 
00272   thePixelsToSkip.clear();
00273   theStDets.theStripsToSkip.clear();
00274 }
00275 
00276 void MeasurementTrackerImpl::updatePixels( const edm::Event& event) const
00277 {
00278   // avoid to update twice from the same event
00279   if (!edm::Service<UpdaterService>()->checkOnce("MeasurementTrackerImpl::updatePixels::"+name_)) return;
00280 
00281   typedef edmNew::DetSet<SiPixelCluster> PixelDetSet;
00282 
00283   bool switchOffPixelsIfEmpty = (!pset_.existsAs<bool>("switchOffPixelsIfEmpty")) ||
00284                                 (pset_.getParameter<bool>("switchOffPixelsIfEmpty"));
00285 
00286   std::vector<uint32_t> rawInactiveDetIds; 
00287   if (!theInactivePixelDetectorLabels.empty()) {
00288     edm::Handle<DetIdCollection> detIds;
00289     for (std::vector<edm::InputTag>::const_iterator itt = theInactivePixelDetectorLabels.begin(), edt = theInactivePixelDetectorLabels.end(); 
00290             itt != edt; ++itt) {
00291       if (event.getByLabel(*itt, detIds)){
00292         rawInactiveDetIds.insert(rawInactiveDetIds.end(), detIds->begin(), detIds->end());
00293       }else{
00294         static bool iFailedAlready=false;
00295         if (!iFailedAlready){
00296           edm::LogError("MissingProduct")<<"I fail to get the list of inactive pixel modules, because of 4.2/4.4 event content change.";
00297           iFailedAlready=true;
00298         }
00299       }
00300     }
00301     if (!rawInactiveDetIds.empty()) std::sort(rawInactiveDetIds.begin(), rawInactiveDetIds.end());
00302   }
00303   // Pixel Clusters
00304   std::string pixelClusterProducer = pset_.getParameter<std::string>("pixelClusterProducer");
00305   if( pixelClusterProducer.empty() ) { //clusters have not been produced
00306     for (std::vector<TkPixelMeasurementDet*>::const_iterator i=thePixelDets.begin();
00307     i!=thePixelDets.end(); i++) {
00308       if (switchOffPixelsIfEmpty) {
00309         (**i).setActiveThisEvent(false);
00310       }else{
00311         (**i).setEmpty();
00312       }
00313     }
00314   }else{  
00315 
00316     edm::Handle<edmNew::DetSetVector<SiPixelCluster> > pixelClusters;
00317     event.getByLabel(pixelClusterProducer, pixelClusters);
00318     const  edmNew::DetSetVector<SiPixelCluster>* pixelCollection = pixelClusters.product();
00319    
00320     if (switchOffPixelsIfEmpty && pixelCollection->empty()) {
00321         for (std::vector<TkPixelMeasurementDet*>::const_iterator i=thePixelDets.begin();
00322              i!=thePixelDets.end(); i++) {
00323               (**i).setActiveThisEvent(false);
00324         }
00325     } else { 
00326 
00327        //std::cout <<"updatePixels "<<pixelCollection->dataSize()<<std::endl;
00328        thePixelsToSkip.resize(pixelCollection->dataSize());
00329        std::fill(thePixelsToSkip.begin(),thePixelsToSkip.end(),false);
00330        
00331        if(selfUpdateSkipClusters_) {
00332           edm::Handle<edm::ContainerMask<edmNew::DetSetVector<SiPixelCluster> > > pixelClusterMask;
00333          //and get the collection of pixel ref to skip
00334          event.getByLabel(pset_.getParameter<edm::InputTag>("skipClusters"),pixelClusterMask);
00335          LogDebug("MeasurementTracker")<<"getting pxl refs to skip";
00336          if (pixelClusterMask.failedToGet())edm::LogError("MeasurementTracker")<<"not getting the pixel clusters to skip";
00337          if (pixelClusterMask->refProd().id()!=pixelClusters.id()){
00338            edm::LogError("ProductIdMismatch")<<"The pixel masking does not point to the proper collection of clusters: "<<pixelClusterMask->refProd().id()<<"!="<<pixelClusters.id();
00339          }
00340          pixelClusterMask->copyMaskTo(thePixelsToSkip);
00341        }
00342           
00343        for (std::vector<TkPixelMeasurementDet*>::const_iterator i=thePixelDets.begin();
00344        i!=thePixelDets.end(); i++) {
00345 
00346          // foreach det get cluster range
00347          unsigned int id = (**i).geomDet().geographicalId().rawId();
00348          if (!rawInactiveDetIds.empty() && std::binary_search(rawInactiveDetIds.begin(), rawInactiveDetIds.end(), id)) {
00349            (**i).setActiveThisEvent(false); continue;
00350          }
00351          //FIXME
00352          //fill the set with what needs to be skipped
00353          edmNew::DetSetVector<SiPixelCluster>::const_iterator it = pixelCollection->find( id );
00354          if ( it != pixelCollection->end() ){            
00355            // push cluster range in det
00356            (**i).update( *it, pixelClusters, id );
00357          } else{
00358            (**i).setEmpty();
00359          }
00360        }
00361     }
00362   }
00363 }
00364 
00365 void MeasurementTrackerImpl::getInactiveStrips(const edm::Event& event,
00366                                                std::vector<uint32_t> & rawInactiveDetIds) const {
00367   if (!theInactiveStripDetectorLabels.empty()) {
00368     edm::Handle<DetIdCollection> detIds;
00369     for (std::vector<edm::InputTag>::const_iterator itt = theInactiveStripDetectorLabels.begin(), edt = theInactiveStripDetectorLabels.end(); 
00370          itt != edt; ++itt) {
00371       event.getByLabel(*itt, detIds);
00372       rawInactiveDetIds.insert(rawInactiveDetIds.end(), detIds->begin(), detIds->end());
00373     }
00374     if (!rawInactiveDetIds.empty()) std::sort(rawInactiveDetIds.begin(), rawInactiveDetIds.end());
00375   }
00376 }
00377 
00378 void MeasurementTrackerImpl::updateStrips( const edm::Event& event) const
00379 {
00380   // avoid to update twice from the same event
00381   if (!edm::Service<UpdaterService>()->checkOnce("MeasurementTrackerImpl::updateStrips::"+name_)) return;
00382 
00383   typedef edmNew::DetSet<SiStripCluster>   StripDetSet;
00384 
00385   std::vector<uint32_t> rawInactiveDetIds;
00386   getInactiveStrips(event,rawInactiveDetIds);
00387 
00388   // Strip Clusters
00389   std::string stripClusterProducer = pset_.getParameter<std::string>("stripClusterProducer");
00390   //first clear all of them
00391   theStDets.setEmpty();
00392 
00393 
00394   if( !stripClusterProducer.compare("") )  return;  //clusters have not been produced
00395 
00396   const int endDet = theStDets.id_.size();
00397  
00398 
00399   // mark as inactive if in rawInactiveDetIds
00400   int i=0;
00401   unsigned int idp=0;
00402   for ( auto id : rawInactiveDetIds) {
00403     if (id==idp) continue; // skip multiple id
00404     idp=id;
00405     i=theStDets.find(id,i);
00406     assert(i!=endDet && id == theStDets.id(i));
00407     theStDets.setActiveThisEvent(i,false);
00408   }
00409 
00410   //=========  actually load cluster =============
00411   if(!theStDets.isRegional()){
00412     edm::Handle<edmNew::DetSetVector<SiStripCluster> > clusterHandle;
00413     event.getByLabel(stripClusterProducer, clusterHandle);
00414     const edmNew::DetSetVector<SiStripCluster>* clusterCollection = clusterHandle.product();
00415     
00416     
00417     if (selfUpdateSkipClusters_){
00418       edm::Handle<edm::ContainerMask<edmNew::DetSetVector<SiStripCluster> > > stripClusterMask;
00419       //and get the collection of pixel ref to skip
00420       LogDebug("MeasurementTracker")<<"getting strp refs to skip";
00421       event.getByLabel(pset_.getParameter<edm::InputTag>("skipClusters"),stripClusterMask);
00422       if (stripClusterMask.failedToGet())  edm::LogError("MeasurementTracker")<<"not getting the strip clusters to skip";
00423       if (stripClusterMask->refProd().id()!=clusterHandle.id()){
00424         edm::LogError("ProductIdMismatch")<<"The strip masking does not point to the proper collection of clusters: "<<stripClusterMask->refProd().id()<<"!="<<clusterHandle.id();
00425       }
00426       stripClusterMask->copyMaskTo(theStDets.theStripsToSkip);
00427     }
00428     
00429     theStDets.handle_ = clusterHandle;
00430     int i=0;
00431     edmNew::DetSetVector<SiStripCluster>::const_iterator it = (*clusterCollection).begin();
00432     edmNew::DetSetVector<SiStripCluster>::const_iterator endColl = (*clusterCollection).end();
00433     // cluster and det and in order (both) and unique so let's use set intersection
00434     for (;it!=endColl; ++it) {
00435       StripDetSet detSet = *it;
00436       unsigned int id = detSet.id();
00437       while ( id != theStDets.id(i)) { // eventually change to lower_bound
00438         ++i;
00439         if (endDet==i) throw "we have a problem!!!!";
00440       }
00441       
00442       // push cluster range in det
00443       if ( theStDets.isActive(i) )
00444         theStDets.update(i,detSet);
00445     }
00446 
00447   }else{   // regional
00448     
00449     //then set the not-empty ones only
00450     edm::Handle<edm::RefGetter<SiStripCluster> > refClusterHandle;
00451     event.getByLabel(stripClusterProducer, refClusterHandle);
00452     
00453     std::string lazyGetter = pset_.getParameter<std::string>("stripLazyGetterProducer");
00454     edm::Handle<edm::LazyGetter<SiStripCluster> > lazyClusterHandle;
00455     event.getByLabel(lazyGetter,lazyClusterHandle);
00456     
00457     if(selfUpdateSkipClusters_){
00458       edm::Handle<edm::ContainerMask<edmNew::DetSetVector<SiStripCluster> > > stripClusterMask;
00459       LogDebug("MeasurementTracker")<<"getting reg strp refs to skip";
00460       event.getByLabel(pset_.getParameter<edm::InputTag>("skipClusters"),stripClusterMask);
00461       if (stripClusterMask.failedToGet())edm::LogError("MeasurementTracker")<<"not getting the strip clusters to skip";
00462       if (stripClusterMask->refProd().id()!=lazyClusterHandle.id()){
00463         edm::LogError("ProductIdMismatch")<<"The strip masking does not point to the proper collection of clusters: "<<stripClusterMask->refProd().id()<<"!="<<lazyClusterHandle.id();
00464       }       
00465       stripClusterMask->copyMaskTo(theStDets.theStripsToSkip);
00466     }
00467     
00468     theStDets.regionalHandle_ =  lazyClusterHandle;
00469     
00470     uint32_t tmpId=0;
00471     vector<SiStripCluster>::const_iterator beginIterator;
00472     edm::RefGetter<SiStripCluster>::const_iterator iregion = refClusterHandle->begin();
00473     for(;iregion!=refClusterHandle->end();++iregion) {
00474       const edm::RegionIndex<SiStripCluster>& region = *iregion;
00475       vector<SiStripCluster>::const_iterator icluster = region.begin();
00476       const vector<SiStripCluster>::const_iterator endIterator = region.end();
00477       tmpId = icluster->geographicalId();
00478       beginIterator = icluster;
00479       
00480       //std::cout << "== tmpId ad inizio loop dentro region: " << tmpId << std::endl;
00481       
00482       for (;icluster!=endIterator;icluster++) {
00483         //std::cout << "===== cluster id,pos " 
00484         //  << icluster->geographicalId() << " , " << icluster->barycenter()
00485         //  << std::endl;
00486         //std::cout << "=====making ref in recHits() " << std::endl;
00487         if( icluster->geographicalId() != tmpId){ 
00488           //std::cout << "geo!=tmpId" << std::endl;
00489           
00490           //cannot we avoid to update the det with detId of itself??  (sure we can!, done!)
00491           theStDets.update(concreteDetUpdatable(tmpId)->index(),beginIterator,icluster);
00492           
00493           tmpId = icluster->geographicalId();
00494           beginIterator = icluster;
00495           if( icluster == (endIterator-1)){
00496             theStDets.update(concreteDetUpdatable(tmpId)->index(),icluster,endIterator);
00497           }   
00498         }else if( icluster == (endIterator-1)){    
00499           theStDets.update(concreteDetUpdatable(tmpId)->index(),beginIterator,endIterator);      
00500         }
00501       }//end loop cluster in one ragion
00502     }
00503   }//end of block for updating with regional clusters 
00504 }
00505 
00506 
00507 TkStripMeasurementDet * MeasurementTrackerImpl::concreteDetUpdatable(DetId id) const {
00508 #ifdef EDM_DEBUG //or similar
00509   const TkStripMeasurementDet* theConcreteDet = 
00510     dynamic_cast<const TkStripMeasurementDet*>(findDet(id));
00511   if(theConcreteDet == 0)
00512     throw MeasurementDetException("failed casting to TkStripMeasurementDet*");      
00513 #endif
00514   // will trigger ondemand unpacking
00515   return const_cast<TkStripMeasurementDet*>(static_cast<const TkStripMeasurementDet*>(idToDet(id)));
00516 }
00517 
00518 
00519 void MeasurementTrackerImpl::initializeStripStatus(const SiStripQuality *quality, int qualityFlags, int qualityDebugFlags) {
00520   edm::ParameterSet cutPset = pset_.getParameter<edm::ParameterSet>("badStripCuts");
00521    theStDets.initializeStripStatus(quality, qualityFlags, qualityDebugFlags, cutPset);
00522 }
00523 
00524 void MeasurementTrackerImpl::initializePixelStatus(const SiPixelQuality *quality, const SiPixelFedCabling *pixelCabling, int qualityFlags, int qualityDebugFlags) {
00525   if ((quality != 0) && (qualityFlags != 0))  {
00526     edm::LogInfo("MeasurementTracker") << "qualityFlags = " << qualityFlags;
00527     unsigned int on = 0, tot = 0, badrocs = 0; 
00528     for (std::vector<TkPixelMeasurementDet*>::const_iterator i=thePixelDets.begin();
00529          i!=thePixelDets.end(); i++) {
00530       uint32_t detid = ((**i).geomDet().geographicalId()).rawId();
00531       if (qualityFlags & BadModules) {
00532           bool isOn = quality->IsModuleUsable(detid);
00533           (*i)->setActive(isOn);
00534           tot++; on += (unsigned int) isOn;
00535           if (qualityDebugFlags & BadModules) {
00536             edm::LogInfo("MeasurementTracker")<< "MeasurementTrackerImpl::initializePixelStatus : detid " << detid << " is " << (isOn ?  "on" : "off");
00537           }
00538        } else {
00539           (*i)->setActive(true);
00540        }
00541        if ((qualityFlags & BadROCs) && (quality->getBadRocs(detid) != 0)) {
00542           std::vector<LocalPoint> badROCs = quality->getBadRocPositions(detid, *theTrackerGeom, pixelCabling);
00543           badrocs += badROCs.size();
00544           (*i)->setBadRocPositions(badROCs);
00545        } else {
00546           (*i)->clearBadRocPositions();  
00547        }
00548     }
00549     if (qualityDebugFlags & BadModules) {
00550         edm::LogInfo("MeasurementTracker PixelModuleStatus") << 
00551             " Total modules: " << tot << ", active " << on <<", inactive " << (tot - on);
00552     }
00553     if (qualityDebugFlags & BadROCs) {
00554         edm::LogInfo("MeasurementTracker PixelROCStatus") << " Total of bad ROCs: " << badrocs ;
00555     }
00556   } else {
00557     for (std::vector<TkPixelMeasurementDet*>::const_iterator i=thePixelDets.begin();
00558          i!=thePixelDets.end(); i++) {
00559       (*i)->setActive(true);          // module ON
00560     }
00561   }
00562 }
00563