CMS 3D CMS Logo

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