CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/RecoLocalTracker/SubCollectionProducers/src/HLTTrackClusterRemover.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/DetSetVectorNew.h"
00021 #include "DataFormats/Common/interface/ContainerMask.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 HLTTrackClusterRemover : public edm::EDProducer {
00035     public:
00036         HLTTrackClusterRemover(const edm::ParameterSet& iConfig) ;
00037         ~HLTTrackClusterRemover() ;
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         bool doStrip_, doPixel_;
00057         edm::InputTag stripClusters_, pixelClusters_;
00058         bool mergeOld_;
00059         edm::InputTag oldRemovalInfo_;
00060 
00061         ParamBlock pblocks_[NumberOfParamBlocks];
00062         void readPSet(const edm::ParameterSet& iConfig, const std::string &name, 
00063                 int id1=-1, int id2=-1, int id3=-1, int id4=-1, int id5=-1, int id6=-1) ;
00064 
00065         std::vector<uint8_t> pixels, strips;                // avoid unneed alloc/dealloc of this
00066         edm::ProductID pixelSourceProdID, stripSourceProdID; // ProdIDs refs must point to (for consistency tests)
00067 
00068         inline void process(const TrackingRecHit *hit, float chi2);
00069         inline void process(const OmniClusterRef & cluRef, uint32_t subdet);
00070 
00071         template<typename T> 
00072         std::auto_ptr<edmNew::DetSetVector<T> >
00073         cleanup(const edmNew::DetSetVector<T> &oldClusters, const std::vector<uint8_t> &isGood, 
00074                     reco::ClusterRemovalInfo::Indices &refs, const reco::ClusterRemovalInfo::Indices *oldRefs) ;
00075 
00076         // Carries in full removal info about a given det from oldRefs
00077         void mergeOld(reco::ClusterRemovalInfo::Indices &refs, const reco::ClusterRemovalInfo::Indices &oldRefs) ;
00078 
00079         bool makeProducts_;
00080         std::vector<bool> collectedRegStrips_;
00081         std::vector<bool> collectedPixels_;        
00082 };
00083 
00084 
00085 using namespace std;
00086 using namespace edm;
00087 using namespace reco;
00088 
00089 void
00090 HLTTrackClusterRemover::readPSet(const edm::ParameterSet& iConfig, const std::string &name, 
00091     int id1, int id2, int id3, int id4, int id5, int id6) 
00092 {
00093     if (iConfig.exists(name)) {
00094         ParamBlock pblock(iConfig.getParameter<ParameterSet>(name));
00095         if (id1 == -1) {
00096             fill(pblocks_, pblocks_+NumberOfParamBlocks, pblock);
00097         } else {
00098             pblocks_[id1] = pblock;
00099             if (id2 != -1) pblocks_[id2] = pblock;
00100             if (id3 != -1) pblocks_[id3] = pblock;
00101             if (id4 != -1) pblocks_[id4] = pblock;
00102             if (id5 != -1) pblocks_[id5] = pblock;
00103             if (id6 != -1) pblocks_[id6] = pblock;
00104         }
00105     }
00106 }
00107 
00108 HLTTrackClusterRemover::HLTTrackClusterRemover(const ParameterSet& iConfig):
00109     trajectories_(iConfig.getParameter<InputTag>("trajectories")),
00110     doStrip_(iConfig.existsAs<bool>("doStrip") ? iConfig.getParameter<bool>("doStrip") : true),
00111     doPixel_(iConfig.existsAs<bool>("doPixel") ? iConfig.getParameter<bool>("doPixel") : true),
00112     stripClusters_(doStrip_ ? iConfig.getParameter<InputTag>("stripClusters") : InputTag("NONE")),
00113     pixelClusters_(doPixel_ ? iConfig.getParameter<InputTag>("pixelClusters") : InputTag("NONE")),
00114     mergeOld_(false),
00115     oldRemovalInfo_(edm::InputTag()),
00116     makeProducts_(true)
00117 {
00118 
00119   if (iConfig.exists("oldClusterRemovalInfo"))
00120     {
00121       oldRemovalInfo_=iConfig.getParameter<InputTag>("oldClusterRemovalInfo");
00122       if (not (oldRemovalInfo_== edm::InputTag())) mergeOld_=true;
00123     }
00124 
00125     fill(pblocks_, pblocks_+NumberOfParamBlocks, ParamBlock());
00126     readPSet(iConfig, "Common",-1);
00127     if (doPixel_) {
00128         readPSet(iConfig, "Pixel" ,0,1);
00129         readPSet(iConfig, "PXB" ,0);
00130         readPSet(iConfig, "PXE" ,1);
00131     }
00132     if (doStrip_) {
00133         readPSet(iConfig, "Strip" ,2,3,4,5);
00134         readPSet(iConfig, "StripInner" ,2,3);
00135         readPSet(iConfig, "StripOuter" ,4,5);
00136         readPSet(iConfig, "TIB" ,2);
00137         readPSet(iConfig, "TID" ,3);
00138         readPSet(iConfig, "TOB" ,4);
00139         readPSet(iConfig, "TEC" ,5);
00140     }
00141 
00142     bool usingCharge = false;
00143     for (size_t i = 0; i < NumberOfParamBlocks; ++i) {
00144         if (!pblocks_[i].isSet_) throw cms::Exception("Configuration Error") << "HLTTrackClusterRemover: Missing configuration for detector with subDetID = " << (i+1);
00145         if (pblocks_[i].usesCharge_ && !usingCharge) {
00146             throw cms::Exception("Configuration Error") << "HLTTrackClusterRemover: Configuration for subDetID = " << (i+1) << " uses cluster charge, which is not enabled.";
00147         }
00148     }
00149 
00150     //produces<edmNew::DetSetVector<SiPixelClusterRefNew> >();
00151     //produces<edmNew::DetSetVector<SiStripRecHit1D::ClusterRegionalRef> >();
00152     
00153     produces< edm::ContainerMask<edmNew::DetSetVector<SiPixelCluster> > >();
00154     produces< edm::ContainerMask<edm::LazyGetter<SiStripCluster> > >();
00155 
00156 }
00157 
00158 
00159 HLTTrackClusterRemover::~HLTTrackClusterRemover()
00160 {
00161 }
00162 
00163 void HLTTrackClusterRemover::mergeOld(ClusterRemovalInfo::Indices &refs,
00164                                             const ClusterRemovalInfo::Indices &oldRefs) 
00165 {
00166         for (size_t i = 0, n = refs.size(); i < n; ++i) {
00167             refs[i] = oldRefs[refs[i]];
00168        }
00169 }
00170 
00171  
00172 template<typename T> 
00173 auto_ptr<edmNew::DetSetVector<T> >
00174 HLTTrackClusterRemover::cleanup(const edmNew::DetSetVector<T> &oldClusters, const std::vector<uint8_t> &isGood, 
00175                              reco::ClusterRemovalInfo::Indices &refs, const reco::ClusterRemovalInfo::Indices *oldRefs){
00176     typedef typename edmNew::DetSetVector<T>             DSV;
00177     typedef typename edmNew::DetSetVector<T>::FastFiller DSF;
00178     typedef typename edmNew::DetSet<T>                   DS;
00179     auto_ptr<DSV> output(new DSV());
00180     output->reserve(oldClusters.size(), oldClusters.dataSize());
00181 
00182     unsigned int countOld=0;
00183     unsigned int countNew=0;
00184 
00185     // cluster removal loop
00186     const T * firstOffset = & oldClusters.data().front();
00187     for (typename DSV::const_iterator itdet = oldClusters.begin(), enddet = oldClusters.end(); itdet != enddet; ++itdet) {
00188         DS oldDS = *itdet;
00189 
00190         if (oldDS.empty()) continue; // skip empty detsets
00191 
00192         uint32_t id = oldDS.detId();
00193         DSF outds(*output, id);
00194 
00195         for (typename DS::const_iterator it = oldDS.begin(), ed = oldDS.end(); it != ed; ++it) {
00196             uint32_t index = ((&*it) - firstOffset);
00197             countOld++;
00198             if (isGood[index]) { 
00199                 outds.push_back(*it);
00200                 countNew++;
00201                 refs.push_back(index); 
00202                 //std::cout << "HLTTrackClusterRemover::cleanup " << typeid(T).name() << " reference " << index << " to " << (refs.size() - 1) << std::endl;
00203             }
00204 
00205         }
00206         if (outds.empty()) outds.abort(); // not write in an empty DSV
00207     }
00208     //    double fraction = countNew  / (double) countOld;
00209     //    std::cout<<"fraction: "<<fraction<<std::endl;
00210     if (oldRefs != 0) mergeOld(refs, *oldRefs);
00211     return output;
00212 }
00213 
00214 void HLTTrackClusterRemover::process(OmniClusterRef const & clusterReg, uint32_t subdet) {
00215    if (clusterReg.id() != stripSourceProdID) throw cms::Exception("Inconsistent Data") <<
00216     "HLTTrackClusterRemover: strip cluster ref from Product ID = " << clusterReg.id() <<
00217     " does not match with source cluster collection (ID = " << stripSourceProdID << ")\n.";
00218   if (collectedRegStrips_.size()<=clusterReg.key()){
00219     edm::LogError("BadCollectionSize")<<collectedRegStrips_.size()<<" is smaller than "<<clusterReg.key();
00220     assert(collectedRegStrips_.size()>clusterReg.key());
00221   }
00222   collectedRegStrips_[clusterReg.key()]=true;
00223 }
00224 
00225 
00226 
00227 void HLTTrackClusterRemover::process(const TrackingRecHit *hit, float chi2) {
00228     DetId detid = hit->geographicalId(); 
00229     uint32_t subdet = detid.subdetId();
00230 
00231     assert ((subdet > 0) && (subdet <= NumberOfParamBlocks));
00232 
00233     // chi2 cut
00234     if (chi2 > pblocks_[subdet-1].maxChi2_) return;
00235 
00236     if ((subdet == PixelSubdetector::PixelBarrel) || (subdet == PixelSubdetector::PixelEndcap)) {
00237       //      std::cout<<"process pxl hit"<<std::endl;
00238         if (!doPixel_) return;
00239         // this is a pixel, and i *know* it is
00240         const SiPixelRecHit *pixelHit = static_cast<const SiPixelRecHit *>(hit);
00241 
00242         SiPixelRecHit::ClusterRef cluster = pixelHit->cluster();
00243 
00244         if (cluster.id() != pixelSourceProdID) throw cms::Exception("Inconsistent Data") << 
00245                 "HLTTrackClusterRemover: pixel cluster ref from Product ID = " << cluster.id() << 
00246                 " does not match with source cluster collection (ID = " << pixelSourceProdID << ")\n.";
00247 
00248         assert(cluster.id() == pixelSourceProdID);
00249 //DBG// cout << "HIT NEW PIXEL DETID = " << detid.rawId() << ", Cluster [ " << cluster.key().first << " / " <<  cluster.key().second << " ] " << endl;
00250 
00251         // if requested, cut on cluster size
00252         if (pblocks_[subdet-1].usesSize_ && (cluster->pixels().size() > pblocks_[subdet-1].maxSize_)) return;
00253 
00254         // mark as used
00255         //pixels[cluster.key()] = false;
00256         assert(collectedPixels_.size() > cluster.key());
00257         collectedPixels_[cluster.key()]=true;
00258     } else { // aka Strip
00259         if (!doStrip_) return;
00260         const type_info &hitType = typeid(*hit);
00261         if (hitType == typeid(SiStripRecHit2D)) {
00262             const SiStripRecHit2D *stripHit = static_cast<const SiStripRecHit2D *>(hit);
00263 //DBG//     cout << "Plain RecHit 2D: " << endl;
00264             process(stripHit->omniClusterRef(),subdet);}
00265         else if (hitType == typeid(SiStripRecHit1D)) {
00266           const SiStripRecHit1D *hit1D = static_cast<const SiStripRecHit1D *>(hit);
00267           process(hit1D->omniClusterRef(),subdet);
00268         } else if (hitType == typeid(SiStripMatchedRecHit2D)) {
00269             const SiStripMatchedRecHit2D *matchHit = static_cast<const SiStripMatchedRecHit2D *>(hit);
00270 //DBG//     cout << "Matched RecHit 2D: " << endl;
00271             process(matchHit->monoClusterRef(),subdet);
00272             process(matchHit->stereoClusterRef(),subdet);
00273         } else if (hitType == typeid(ProjectedSiStripRecHit2D)) {
00274             const ProjectedSiStripRecHit2D *projHit = static_cast<const ProjectedSiStripRecHit2D *>(hit);
00275 //DBG//     cout << "Projected RecHit 2D: " << endl;
00276             process(projHit->originalHit().omniClusterRef(),subdet);
00277         } else throw cms::Exception("NOT IMPLEMENTED") << "Don't know how to handle " << hitType.name() << " on detid " << detid.rawId() << "\n";
00278     }
00279 }
00280 
00281 /*   Schematic picture of n-th step Iterative removal
00282  *   (that os removing clusters after n-th step tracking)
00283  *   clusters:    [ C1 ] -> [ C2 ] -> ... -> [ Cn ] -> [ Cn + 1 ]
00284  *                   ^                          ^           ^--- OUTPUT "new" ID
00285  *                   |-- before any removal     |----- Source clusters                                                                             
00286  *                   |-- OUTPUT "old" ID        |----- Hits in Traj. point here                       
00287  *                   |                          \----- Old ClusterRemovalInfo "new" ID                 
00288  *                   \-- Old ClusterRemovalInfo "old" ID                                                  
00289  */
00290 
00291 
00292 void
00293 HLTTrackClusterRemover::produce(Event& iEvent, const EventSetup& iSetup)
00294 {
00295     ProductID pixelOldProdID, stripOldProdID;
00296 
00297     Handle<edmNew::DetSetVector<SiPixelCluster> > pixelClusters;
00298     if (doPixel_) {
00299         iEvent.getByLabel(pixelClusters_, pixelClusters);
00300         pixelSourceProdID = pixelClusters.id();
00301     }
00302 
00303     Handle<edm::LazyGetter<SiStripCluster> > stripClusters;
00304     if (doStrip_) {
00305         iEvent.getByLabel(stripClusters_, stripClusters);
00306         stripSourceProdID = stripClusters.id();
00307     }
00308 
00309     //Handle<TrajTrackAssociationCollection> trajectories; 
00310     Handle<vector<Trajectory> > trajectories; 
00311     iEvent.getByLabel(trajectories_, trajectories);
00312 
00313     typedef edm::ContainerMask<edmNew::DetSetVector<SiPixelCluster> > PixelMaskContainer;
00314     typedef edm::ContainerMask<edm::LazyGetter<SiStripCluster> > StripMaskContainer;
00315     if(mergeOld_) {
00316       edm::Handle<PixelMaskContainer> oldPxlMask;
00317       edm::Handle<StripMaskContainer> oldStrMask;
00318       iEvent.getByLabel(oldRemovalInfo_,oldPxlMask);
00319       iEvent.getByLabel(oldRemovalInfo_,oldStrMask);
00320       LogDebug("TrackClusterRemover")<<"to merge in, "<<oldStrMask->size()<<" strp and "<<oldPxlMask->size()<<" pxl";
00321       oldStrMask->copyMaskTo(collectedRegStrips_);
00322       oldPxlMask->copyMaskTo(collectedPixels_);
00323       collectedRegStrips_.resize(stripClusters->size());
00324     }else {
00325       collectedRegStrips_.resize(stripClusters->size()); fill(collectedRegStrips_.begin(), collectedRegStrips_.end(), false);
00326       collectedPixels_.resize(pixelClusters->dataSize()); fill(collectedPixels_.begin(), collectedPixels_.end(), false);
00327     } 
00328 
00329 
00330     //for (TrajTrackAssociationCollection::const_iterator it = trajectories->begin(), ed = trajectories->end(); it != ed; ++it)  {
00331     //    const Trajectory &tj = * it->key;
00332     for (vector<Trajectory>::const_iterator it = trajectories->begin(), ed = trajectories->end(); it != ed; ++it)  {
00333         const Trajectory &tj = * it;
00334         const vector<TrajectoryMeasurement> &tms = tj.measurements();
00335         vector<TrajectoryMeasurement>::const_iterator itm, endtm;
00336         for (itm = tms.begin(), endtm = tms.end(); itm != endtm; ++itm) {
00337             const TrackingRecHit *hit = itm->recHit()->hit();
00338             if (!hit->isValid()) continue; 
00339             //      std::cout<<"process hit"<<std::endl;
00340             process( hit, itm->estimate() );
00341         }
00342     }
00343     
00344     std::auto_ptr<StripMaskContainer> removedStripClusterMask(
00345        new StripMaskContainer(edm::RefProd<edm::LazyGetter<SiStripCluster> >(stripClusters),collectedRegStrips_));
00346     LogDebug("TrackClusterRemover")<<"total strip to skip: "<<std::count(collectedRegStrips_.begin(),collectedRegStrips_.end(),true);
00347     iEvent.put( removedStripClusterMask );
00348 
00349     std::auto_ptr<PixelMaskContainer> removedPixelClusterMask(
00350        new PixelMaskContainer(edm::RefProd<edmNew::DetSetVector<SiPixelCluster> >(pixelClusters),collectedPixels_));      
00351     LogDebug("TrackClusterRemover")<<"total pxl to skip: "<<std::count(collectedPixels_.begin(),collectedPixels_.end(),true);
00352     iEvent.put( removedPixelClusterMask );
00353     
00354 }
00355 
00356 #include "FWCore/PluginManager/interface/ModuleDef.h"
00357 #include "FWCore/Framework/interface/MakerMacros.h"
00358 DEFINE_FWK_MODULE(HLTTrackClusterRemover);