CMS 3D CMS Logo

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