CMS 3D CMS Logo

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