CMS 3D CMS Logo

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