CMS 3D CMS Logo

CMSSW_4_4_3_patch1/src/RecoLocalTracker/SubCollectionProducers/src/TrackClusterRemover.cc

Go to the documentation of this file.
00001 #include "FWCore/Framework/interface/Frameworkfwd.h"
00002 #include "FWCore/Framework/interface/EDProducer.h"
00003 #include "FWCore/Framework/interface/Event.h"
00004 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00005 #include "FWCore/Utilities/interface/InputTag.h"
00006 
00007 #include "DataFormats/SiStripCluster/interface/SiStripCluster.h"
00008 #include "DataFormats/SiPixelCluster/interface/SiPixelCluster.h"
00009 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2D.h"
00010 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit1D.h"
00011 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2D.h"
00012 #include "DataFormats/TrackerRecHit2D/interface/ProjectedSiStripRecHit2D.h"
00013 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
00014 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHitCollection.h"
00015 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2DCollection.h"
00016 #include "DataFormats/Common/interface/Handle.h"
00017 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
00018 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
00019 #include "DataFormats/Common/interface/DetSetVector.h"
00020 #include "DataFormats/Common/interface/ValueMap.h"
00021 #include "DataFormats/Common/interface/DetSetVectorNew.h"
00022 #include "DataFormats/Provenance/interface/ProductID.h"
00023 
00024 #include "DataFormats/TrackReco/interface/Track.h"
00025 #include "DataFormats/TrackerRecHit2D/interface/ClusterRemovalInfo.h"
00026 
00027 #include "TrackingTools/PatternTools/interface/Trajectory.h"
00028 #include "TrackingTools/PatternTools/interface/TrajTrackAssociation.h"
00029 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00030 //
00031 // class decleration
00032 //
00033 
00034 class TrackClusterRemover : public edm::EDProducer {
00035     public:
00036         TrackClusterRemover(const edm::ParameterSet& iConfig) ;
00037         ~TrackClusterRemover() ;
00038         void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) ;
00039     private:
00040         struct ParamBlock {
00041             ParamBlock() : isSet_(false), usesCharge_(false) {}
00042             ParamBlock(const edm::ParameterSet& iConfig) :
00043                 isSet_(true), 
00044                 usesCharge_(iConfig.exists("maxCharge")),
00045                 usesSize_(iConfig.exists("maxSize")),
00046                 maxChi2_(iConfig.getParameter<double>("maxChi2")),
00047                 maxCharge_(usesCharge_ ? iConfig.getParameter<double>("maxCharge") : 0), 
00048                 maxSize_(usesSize_ ? iConfig.getParameter<uint32_t>("maxSize") : 0) { }
00049             bool  isSet_, usesCharge_, usesSize_;
00050             float maxChi2_, maxCharge_;
00051             size_t maxSize_;
00052         };
00053         static const unsigned int NumberOfParamBlocks = 6;
00054 
00055         edm::InputTag trajectories_;
00056         std::vector<edm::InputTag> overrideTrkQuals_;
00057         bool doStrip_, doPixel_;
00058         edm::InputTag stripClusters_, pixelClusters_;
00059         bool mergeOld_;
00060         edm::InputTag oldRemovalInfo_;
00061 
00062         ParamBlock pblocks_[NumberOfParamBlocks];
00063         void readPSet(const edm::ParameterSet& iConfig, const std::string &name, 
00064                 int id1=-1, int id2=-1, int id3=-1, int id4=-1, int id5=-1, int id6=-1) ;
00065 
00066         std::vector<uint8_t> pixels, strips;                // avoid unneed alloc/dealloc of this
00067         edm::ProductID pixelSourceProdID, stripSourceProdID; // ProdIDs refs must point to (for consistency tests)
00068 
00069         inline void process(const TrackingRecHit *hit, float chi2);
00070         inline void process(const SiStripRecHit2D *hit2D, uint32_t subdet);
00071         inline void process(const SiStripRecHit1D *hit1D, uint32_t subdet);
00072 
00073         template<typename T> 
00074         std::auto_ptr<edmNew::DetSetVector<T> >
00075         cleanup(const edmNew::DetSetVector<T> &oldClusters, const std::vector<uint8_t> &isGood, 
00076                     reco::ClusterRemovalInfo::Indices &refs, const reco::ClusterRemovalInfo::Indices *oldRefs) ;
00077 
00078         // Carries in full removal info about a given det from oldRefs
00079         void mergeOld(reco::ClusterRemovalInfo::Indices &refs, const reco::ClusterRemovalInfo::Indices &oldRefs) ;
00080 
00081   bool clusterWasteSolution_;
00082   bool filterTracks_;
00083   reco::TrackBase::TrackQuality trackQuality_;
00084   std::map<uint32_t, std::set< SiStripRecHit1D::ClusterRef > > collectedStrip;
00085   std::map<uint32_t, std::set< SiPixelRecHit::ClusterRef  > > collectedPixel;
00086 };
00087 
00088 
00089 using namespace std;
00090 using namespace edm;
00091 using namespace reco;
00092 
00093 void
00094 TrackClusterRemover::readPSet(const edm::ParameterSet& iConfig, const std::string &name, 
00095     int id1, int id2, int id3, int id4, int id5, int id6) 
00096 {
00097     if (iConfig.exists(name)) {
00098         ParamBlock pblock(iConfig.getParameter<ParameterSet>(name));
00099         if (id1 == -1) {
00100             fill(pblocks_, pblocks_+NumberOfParamBlocks, pblock);
00101         } else {
00102             pblocks_[id1] = pblock;
00103             if (id2 != -1) pblocks_[id2] = pblock;
00104             if (id3 != -1) pblocks_[id3] = pblock;
00105             if (id4 != -1) pblocks_[id4] = pblock;
00106             if (id5 != -1) pblocks_[id5] = pblock;
00107             if (id6 != -1) pblocks_[id6] = pblock;
00108         }
00109     }
00110 }
00111 
00112 TrackClusterRemover::TrackClusterRemover(const ParameterSet& iConfig):
00113     trajectories_(iConfig.getParameter<InputTag>("trajectories")),
00114     doStrip_(iConfig.existsAs<bool>("doStrip") ? iConfig.getParameter<bool>("doStrip") : true),
00115     doPixel_(iConfig.existsAs<bool>("doPixel") ? iConfig.getParameter<bool>("doPixel") : true),
00116     stripClusters_(doStrip_ ? iConfig.getParameter<InputTag>("stripClusters") : InputTag("NONE")),
00117     pixelClusters_(doPixel_ ? iConfig.getParameter<InputTag>("pixelClusters") : InputTag("NONE")),
00118     mergeOld_(iConfig.exists("oldClusterRemovalInfo")),
00119     oldRemovalInfo_(mergeOld_ ? iConfig.getParameter<InputTag>("oldClusterRemovalInfo") : InputTag("NONE")),
00120     clusterWasteSolution_(true)
00121 {
00122   if (iConfig.exists("overrideTrkQuals"))
00123     overrideTrkQuals_.push_back(iConfig.getParameter<edm::InputTag>("overrideTrkQuals"));
00124   if (iConfig.exists("clusterLessSolution"))
00125     clusterWasteSolution_=!iConfig.getParameter<bool>("clusterLessSolution");
00126   if (doPixel_ && clusterWasteSolution_) produces< edmNew::DetSetVector<SiPixelCluster> >();
00127   if (doStrip_ && clusterWasteSolution_) produces< edmNew::DetSetVector<SiStripCluster> >();
00128   if (clusterWasteSolution_) produces< ClusterRemovalInfo >();
00129 
00130     fill(pblocks_, pblocks_+NumberOfParamBlocks, ParamBlock());
00131     readPSet(iConfig, "Common",-1);
00132     if (doPixel_) {
00133         readPSet(iConfig, "Pixel" ,0,1);
00134         readPSet(iConfig, "PXB" ,0);
00135         readPSet(iConfig, "PXE" ,1);
00136     }
00137     if (doStrip_) {
00138         readPSet(iConfig, "Strip" ,2,3,4,5);
00139         readPSet(iConfig, "StripInner" ,2,3);
00140         readPSet(iConfig, "StripOuter" ,4,5);
00141         readPSet(iConfig, "TIB" ,2);
00142         readPSet(iConfig, "TID" ,3);
00143         readPSet(iConfig, "TOB" ,4);
00144         readPSet(iConfig, "TEC" ,5);
00145     }
00146 
00147     bool usingCharge = false;
00148     for (size_t i = 0; i < NumberOfParamBlocks; ++i) {
00149         if (!pblocks_[i].isSet_) throw cms::Exception("Configuration Error") << "TrackClusterRemover: Missing configuration for detector with subDetID = " << (i+1);
00150         if (pblocks_[i].usesCharge_ && !usingCharge) {
00151             throw cms::Exception("Configuration Error") << "TrackClusterRemover: Configuration for subDetID = " << (i+1) << " uses cluster charge, which is not enabled.";
00152         }
00153     }
00154 
00155     if (!clusterWasteSolution_){
00156       produces<edmNew::DetSetVector<SiPixelClusterRefNew> >();
00157       produces<edmNew::DetSetVector<SiStripRecHit1D::ClusterRef> >();
00158     }
00159     trackQuality_=reco::TrackBase::undefQuality;
00160     filterTracks_=false;
00161     if (iConfig.exists("TrackQuality")){
00162       filterTracks_=true;
00163       trackQuality_=reco::TrackBase::qualityByName(iConfig.getParameter<std::string>("TrackQuality"));
00164     }
00165 
00166 }
00167 
00168 
00169 TrackClusterRemover::~TrackClusterRemover()
00170 {
00171 }
00172 
00173 void TrackClusterRemover::mergeOld(ClusterRemovalInfo::Indices &refs,
00174                                             const ClusterRemovalInfo::Indices &oldRefs) 
00175 {
00176         for (size_t i = 0, n = refs.size(); i < n; ++i) {
00177             refs[i] = oldRefs[refs[i]];
00178        }
00179 }
00180 
00181  
00182 template<typename T> 
00183 auto_ptr<edmNew::DetSetVector<T> >
00184 TrackClusterRemover::cleanup(const edmNew::DetSetVector<T> &oldClusters, const std::vector<uint8_t> &isGood, 
00185                              reco::ClusterRemovalInfo::Indices &refs, const reco::ClusterRemovalInfo::Indices *oldRefs){
00186     typedef typename edmNew::DetSetVector<T>             DSV;
00187     typedef typename edmNew::DetSetVector<T>::FastFiller DSF;
00188     typedef typename edmNew::DetSet<T>                   DS;
00189     auto_ptr<DSV> output(new DSV());
00190     output->reserve(oldClusters.size(), oldClusters.dataSize());
00191 
00192     unsigned int countOld=0;
00193     unsigned int countNew=0;
00194 
00195     // cluster removal loop
00196     const T * firstOffset = & oldClusters.data().front();
00197     for (typename DSV::const_iterator itdet = oldClusters.begin(), enddet = oldClusters.end(); itdet != enddet; ++itdet) {
00198         DS oldDS = *itdet;
00199 
00200         if (oldDS.empty()) continue; // skip empty detsets
00201 
00202         uint32_t id = oldDS.detId();
00203         DSF outds(*output, id);
00204 
00205         for (typename DS::const_iterator it = oldDS.begin(), ed = oldDS.end(); it != ed; ++it) {
00206             uint32_t index = ((&*it) - firstOffset);
00207             countOld++;
00208             if (isGood[index]) { 
00209                 outds.push_back(*it);
00210                 countNew++;
00211                 refs.push_back(index); 
00212                 //std::cout << "TrackClusterRemover::cleanup " << typeid(T).name() << " reference " << index << " to " << (refs.size() - 1) << std::endl;
00213             }
00214 
00215         }
00216         if (outds.empty()) outds.abort(); // not write in an empty DSV
00217     }
00218 
00219     if (oldRefs != 0) mergeOld(refs, *oldRefs);
00220     return output;
00221 }
00222 
00223 
00224 void TrackClusterRemover::process(const SiStripRecHit2D *hit, uint32_t subdet) {
00225   SiStripRecHit2D::ClusterRef cluster = hit->cluster();
00226   if (cluster.id() != stripSourceProdID) throw cms::Exception("Inconsistent Data") <<
00227     "TrackClusterRemover: strip cluster ref from Product ID = " << cluster.id() <<
00228     " does not match with source cluster collection (ID = " << stripSourceProdID << ")\n.";
00229   
00230   assert(cluster.id() == stripSourceProdID);
00231   if (pblocks_[subdet-1].usesSize_ && (cluster->amplitudes().size() > pblocks_[subdet-1].maxSize_)) return;
00232 
00233   strips[cluster.key()] = false;  
00234   if (!clusterWasteSolution_) collectedStrip[hit->geographicalId()].insert(cluster);
00235 }
00236 
00237 void TrackClusterRemover::process(const SiStripRecHit1D *hit, uint32_t subdet) {
00238     SiStripRecHit1D::ClusterRef cluster = hit->cluster();
00239     if (cluster.id() != stripSourceProdID) throw cms::Exception("Inconsistent Data") << 
00240             "TrackClusterRemover: strip cluster ref from Product ID = " << cluster.id() << 
00241             " does not match with source cluster collection (ID = " << stripSourceProdID << ")\n.";
00242 
00243     assert(cluster.id() == stripSourceProdID);
00244     if (pblocks_[subdet-1].usesSize_ && (cluster->amplitudes().size() > pblocks_[subdet-1].maxSize_)) return;
00245 
00246 //DBG// cout << "Individual HIT " << cluster.key().first << ", INDEX = " << cluster.key().second << endl;
00247     strips[cluster.key()] = false;
00248     if (!clusterWasteSolution_) collectedStrip[hit->geographicalId()].insert(cluster);
00249 }
00250 
00251 
00252 
00253 
00254 void TrackClusterRemover::process(const TrackingRecHit *hit, float chi2) {
00255     DetId detid = hit->geographicalId(); 
00256     uint32_t subdet = detid.subdetId();
00257 
00258     assert ((subdet > 0) && (subdet <= NumberOfParamBlocks));
00259 
00260     // chi2 cut
00261     if (chi2 > pblocks_[subdet-1].maxChi2_) return;
00262 
00263     if ((subdet == PixelSubdetector::PixelBarrel) || (subdet == PixelSubdetector::PixelEndcap)) {
00264         if (!doPixel_) return;
00265         // this is a pixel, and i *know* it is
00266         const SiPixelRecHit *pixelHit = static_cast<const SiPixelRecHit *>(hit);
00267 
00268         SiPixelRecHit::ClusterRef cluster = pixelHit->cluster();
00269 
00270         if (cluster.id() != pixelSourceProdID) throw cms::Exception("Inconsistent Data") << 
00271                 "TrackClusterRemover: pixel cluster ref from Product ID = " << cluster.id() << 
00272                 " does not match with source cluster collection (ID = " << pixelSourceProdID << ")\n.";
00273 
00274         assert(cluster.id() == pixelSourceProdID);
00275 //DBG// cout << "HIT NEW PIXEL DETID = " << detid.rawId() << ", Cluster [ " << cluster.key().first << " / " <<  cluster.key().second << " ] " << endl;
00276 
00277         // if requested, cut on cluster size
00278         if (pblocks_[subdet-1].usesSize_ && (cluster->pixels().size() > pblocks_[subdet-1].maxSize_)) return;
00279 
00280         // mark as used
00281         pixels[cluster.key()] = false;
00282         
00283         if(!clusterWasteSolution_) collectedPixel[detid.rawId()].insert(cluster);
00284     } else { // aka Strip
00285         if (!doStrip_) return;
00286         const type_info &hitType = typeid(*hit);
00287         if (hitType == typeid(SiStripRecHit2D)) {
00288             const SiStripRecHit2D *stripHit = static_cast<const SiStripRecHit2D *>(hit);
00289 //DBG//     cout << "Plain RecHit 2D: " << endl;
00290             process(stripHit,subdet);}
00291         else if (hitType == typeid(SiStripRecHit1D)) {
00292           const SiStripRecHit1D *hit1D = static_cast<const SiStripRecHit1D *>(hit);
00293           process(hit1D,subdet);
00294         } else if (hitType == typeid(SiStripMatchedRecHit2D)) {
00295             const SiStripMatchedRecHit2D *matchHit = static_cast<const SiStripMatchedRecHit2D *>(hit);
00296 //DBG//     cout << "Matched RecHit 2D: " << endl;
00297             process(matchHit->monoHit(),subdet);
00298             process(matchHit->stereoHit(),subdet);
00299         } else if (hitType == typeid(ProjectedSiStripRecHit2D)) {
00300             const ProjectedSiStripRecHit2D *projHit = static_cast<const ProjectedSiStripRecHit2D *>(hit);
00301 //DBG//     cout << "Projected RecHit 2D: " << endl;
00302             process(&projHit->originalHit(),subdet);
00303         } else throw cms::Exception("NOT IMPLEMENTED") << "Don't know how to handle " << hitType.name() << " on detid " << detid.rawId() << "\n";
00304     }
00305 }
00306 
00307 /*   Schematic picture of n-th step Iterative removal
00308  *   (that os removing clusters after n-th step tracking)
00309  *   clusters:    [ C1 ] -> [ C2 ] -> ... -> [ Cn ] -> [ Cn + 1 ]
00310  *                   ^                          ^           ^--- OUTPUT "new" ID
00311  *                   |-- before any removal     |----- Source clusters                                                                             
00312  *                   |-- OUTPUT "old" ID        |----- Hits in Traj. point here                       
00313  *                   |                          \----- Old ClusterRemovalInfo "new" ID                 
00314  *                   \-- Old ClusterRemovalInfo "old" ID                                                  
00315  */
00316 
00317 
00318 void
00319 TrackClusterRemover::produce(Event& iEvent, const EventSetup& iSetup)
00320 {
00321     ProductID pixelOldProdID, stripOldProdID;
00322 
00323     Handle<edmNew::DetSetVector<SiPixelCluster> > pixelClusters;
00324     if (doPixel_) {
00325         iEvent.getByLabel(pixelClusters_, pixelClusters);
00326         pixelSourceProdID = pixelClusters.id();
00327     }
00328 //DBG// std::cout << "TrackClusterRemover: Read pixel " << pixelClusters_.encode() << " = ID " << pixelSourceProdID << std::endl;
00329 
00330     Handle<edmNew::DetSetVector<SiStripCluster> > stripClusters;
00331     if (doStrip_) {
00332         iEvent.getByLabel(stripClusters_, stripClusters);
00333         stripSourceProdID = stripClusters.id();
00334     }
00335 //DBG// std::cout << "TrackClusterRemover: Read strip " << stripClusters_.encode() << " = ID " << stripSourceProdID << std::endl;
00336 
00337     auto_ptr<ClusterRemovalInfo> cri;
00338     if (clusterWasteSolution_){
00339       if (doStrip_ && doPixel_) cri.reset(new ClusterRemovalInfo(pixelClusters, stripClusters));
00340       else if (doStrip_) cri.reset(new ClusterRemovalInfo(stripClusters));
00341       else if (doPixel_) cri.reset(new ClusterRemovalInfo(pixelClusters));
00342     }
00343 
00344     Handle<ClusterRemovalInfo> oldRemovalInfo;
00345     if (mergeOld_ && clusterWasteSolution_) { 
00346         iEvent.getByLabel(oldRemovalInfo_, oldRemovalInfo); 
00347         // Check ProductIDs
00348         if ( (oldRemovalInfo->stripNewRefProd().id() == stripClusters.id()) &&
00349              (oldRemovalInfo->pixelNewRefProd().id() == pixelClusters.id()) ) {
00350 
00351             cri->getOldClustersFrom(*oldRemovalInfo);
00352 
00353             pixelOldProdID = oldRemovalInfo->pixelRefProd().id();
00354             stripOldProdID = oldRemovalInfo->stripRefProd().id();
00355 
00356         } else {
00357             throw cms::Exception("Inconsistent Data") << "TrackClusterRemover: " <<
00358                 "Input collection product IDs are [pixel: " << pixelClusters.id() << ", strip: " << stripClusters.id() << "] \n" <<
00359                 "\t but the *old* ClusterRemovalInfo " << oldRemovalInfo_.encode() << " refers as 'new product ids' to " <<
00360                     "[pixel: " << oldRemovalInfo->pixelNewRefProd().id() << ", strip: " << oldRemovalInfo->stripNewRefProd().id() << "]\n" << 
00361                 "NOTA BENE: when running TrackClusterRemover with an old ClusterRemovalInfo the hits in the trajectory MUST be already re-keyed.\n";
00362         }
00363     } else { // then Old == Source
00364         pixelOldProdID = pixelSourceProdID;
00365         stripOldProdID = stripSourceProdID;
00366     }
00367 
00368     Handle<TrajTrackAssociationCollection> trajectories_totrack; 
00369     iEvent.getByLabel(trajectories_,trajectories_totrack);
00370 
00371     std::vector<Handle<edm::ValueMap<int> > > quals;
00372     if ( overrideTrkQuals_.size() > 0) {
00373       quals.resize(1);
00374       iEvent.getByLabel(overrideTrkQuals_[0],quals[0]);
00375     }
00376 
00377     if (doStrip_) {
00378       strips.resize(stripClusters->dataSize()); fill(strips.begin(), strips.end(), true);
00379     }
00380     if (doPixel_) {
00381       pixels.resize(pixelClusters->dataSize()); fill(pixels.begin(), pixels.end(), true);
00382     }
00383 
00384     TrajTrackAssociationCollection::const_iterator asst=trajectories_totrack->begin();
00385    
00386     for ( ; asst!=trajectories_totrack->end();++asst){
00387       const Track & track = *(asst->val);
00388       if (filterTracks_) {
00389         bool goodTk = true;
00390         if ( quals.size()!=0) {
00391           int qual=(*(quals[0]))[asst->val];
00392           if ( qual < 0 ) {goodTk=false;}
00393           //note that this does not work for some trackquals (goodIterative  or undefQuality)
00394           else
00395             goodTk = ( qual & (1<<trackQuality_))>>trackQuality_;
00396         }
00397         else
00398           goodTk=(track.quality(trackQuality_));
00399         if ( !goodTk) continue;
00400       }
00401       const Trajectory &tj = *(asst->key);
00402       const vector<TrajectoryMeasurement> &tms = tj.measurements();
00403       vector<TrajectoryMeasurement>::const_iterator itm, endtm;
00404       for (itm = tms.begin(), endtm = tms.end(); itm != endtm; ++itm) {
00405         const TrackingRecHit *hit = itm->recHit()->hit();
00406         if (!hit->isValid()) continue; 
00407         process( hit, itm->estimate() );
00408       }
00409     }
00410 
00411     
00412     if (doPixel_ && clusterWasteSolution_) {
00413         auto_ptr<edmNew::DetSetVector<SiPixelCluster> > newPixelClusters = cleanup(*pixelClusters, pixels, 
00414                     cri->pixelIndices(), mergeOld_ ? &oldRemovalInfo->pixelIndices() : 0);
00415         OrphanHandle<edmNew::DetSetVector<SiPixelCluster> > newPixels = iEvent.put(newPixelClusters); 
00416 //DBG// std::cout << "TrackClusterRemover: Wrote pixel " << newPixels.id() << " from " << pixelSourceProdID << std::endl;
00417         cri->setNewPixelClusters(newPixels);
00418     }
00419     if (doStrip_ && clusterWasteSolution_) {
00420         auto_ptr<edmNew::DetSetVector<SiStripCluster> > newStripClusters = cleanup(*stripClusters, strips, 
00421                     cri->stripIndices(), mergeOld_ ? &oldRemovalInfo->stripIndices() : 0);
00422         OrphanHandle<edmNew::DetSetVector<SiStripCluster> > newStrips = iEvent.put(newStripClusters); 
00423 //DBG// std::cout << "TrackClusterRemover: Wrote strip " << newStrips.id() << " from " << stripSourceProdID << std::endl;
00424         cri->setNewStripClusters(newStrips);
00425     }
00426 
00427     
00428     if (clusterWasteSolution_) {
00429       //      double fraction_pxl= cri->pixelIndices().size() / (double) pixels.size();
00430       //      double fraction_strp= cri->stripIndices().size() / (double) strips.size();
00431       //      edm::LogWarning("TrackClusterRemover")<<" fraction: " << fraction_pxl <<" "<<fraction_strp;
00432       iEvent.put(cri);
00433     }
00434 
00435     pixels.clear(); strips.clear(); 
00436 
00437     if (!clusterWasteSolution_){
00438       auto_ptr<edmNew::DetSetVector<SiPixelClusterRefNew> > removedPixelClsuterRefs(new edmNew::DetSetVector<SiPixelClusterRefNew>());
00439       auto_ptr<edmNew::DetSetVector<SiStripRecHit1D::ClusterRef> > removedStripClsuterRefs(new edmNew::DetSetVector<SiStripRecHit1D::ClusterRef>());
00440       
00441       edm::Handle<edmNew::DetSetVector<SiPixelClusterRefNew> > oldPxlRef;
00442       edm::Handle<edmNew::DetSetVector<SiStripRecHit1D::ClusterRef> > oldStrRef;
00443       if (mergeOld_){
00444         iEvent.getByLabel(oldRemovalInfo_,oldPxlRef);
00445         iEvent.getByLabel(oldRemovalInfo_,oldStrRef);
00446         LogDebug("TrackClusterRemover")<<"to merge in, "<<oldStrRef->size()<<" strp and "<<oldPxlRef->size()<<" pxl";
00447       }      
00448 
00449       if (mergeOld_){
00450         for (edmNew::DetSetVector<SiStripRecHit1D::ClusterRef>::const_iterator itOld=oldStrRef->begin();itOld!=oldStrRef->end();++itOld){
00451           uint32_t id = itOld->detId();
00452           collectedStrip[id].insert(itOld->begin(),itOld->end());
00453         }
00454       }
00455       for (std::map<uint32_t, std::set< SiStripRecHit1D::ClusterRef > >::iterator itskiped= collectedStrip.begin();
00456            itskiped!=collectedStrip.end();++itskiped){
00457         edmNew::DetSetVector<SiStripRecHit1D::ClusterRef>::FastFiller fill(*removedStripClsuterRefs, itskiped->first);
00458         for (std::set< SiStripRecHit1D::ClusterRef >::iterator topush = itskiped->second.begin();
00459              topush!=itskiped->second.end();++topush){
00460           fill.push_back(*topush);
00461           LogDebug("TrackClusterRemover")<<"registering strp ref to be skip on: "<<itskiped->first<<" key: "<<topush->key();
00462         }
00463         if (fill.empty()) fill.abort();
00464         itskiped->second.clear();
00465       }
00466       LogDebug("TrackClusterRemover")<<"total strip to skip: "<<removedStripClsuterRefs->size();
00467       iEvent.put( removedStripClsuterRefs );
00468 
00469       if (mergeOld_){
00470         for (edmNew::DetSetVector<SiPixelClusterRefNew>::const_iterator itOld=oldPxlRef->begin();itOld!=oldPxlRef->end();++itOld){
00471           uint32_t id = itOld->detId();
00472           collectedPixel[id].insert(itOld->begin(),itOld->end());
00473         }
00474       }
00475 
00476       for (std::map<uint32_t, std::set< SiPixelRecHit::ClusterRef  > >::iterator itskiped= collectedPixel.begin();
00477            itskiped!=collectedPixel.end();++itskiped){  
00478         edmNew::DetSetVector<SiPixelClusterRefNew>::FastFiller fill(*removedPixelClsuterRefs, itskiped->first);
00479         for (std::set< SiPixelRecHit::ClusterRef  >::iterator topush = itskiped->second.begin(); 
00480              topush!=itskiped->second.end();++topush){ 
00481           fill.push_back(*topush); 
00482           LogDebug("TrackClusterRemover")<<"registering pxk ref to be skip on: "<<itskiped->first<<" key: "<<topush->key();
00483         }
00484         if (fill.empty()) fill.abort();   
00485         itskiped->second.clear();
00486       }
00487       LogDebug("TrackClusterRemover")<<"total pxl to skip: "<<removedPixelClsuterRefs->size();
00488       iEvent.put( removedPixelClsuterRefs );
00489       
00490     }
00491     
00492 
00493 }
00494 
00495 #include "FWCore/PluginManager/interface/ModuleDef.h"
00496 #include "FWCore/Framework/interface/MakerMacros.h"
00497 DEFINE_FWK_MODULE(TrackClusterRemover);