CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch13/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/DetSetVectorNew.h"
00021 #include "DataFormats/Provenance/interface/ProductID.h"
00022 
00023 #include "DataFormats/TrackReco/interface/Track.h"
00024 #include "DataFormats/TrackerRecHit2D/interface/ClusterRemovalInfo.h"
00025 
00026 #include "TrackingTools/PatternTools/interface/Trajectory.h"
00027 #include "TrackingTools/PatternTools/interface/TrajTrackAssociation.h"
00028 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00029 //
00030 // class decleration
00031 //
00032 
00033 class TrackClusterRemover : public edm::EDProducer {
00034     public:
00035         TrackClusterRemover(const edm::ParameterSet& iConfig) ;
00036         ~TrackClusterRemover() ;
00037         void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) ;
00038     private:
00039         struct ParamBlock {
00040             ParamBlock() : isSet_(false), usesCharge_(false) {}
00041             ParamBlock(const edm::ParameterSet& iConfig) :
00042                 isSet_(true), 
00043                 usesCharge_(iConfig.exists("maxCharge")),
00044                 usesSize_(iConfig.exists("maxSize")),
00045                 maxChi2_(iConfig.getParameter<double>("maxChi2")),
00046                 maxCharge_(usesCharge_ ? iConfig.getParameter<double>("maxCharge") : 0), 
00047                 maxSize_(usesSize_ ? iConfig.getParameter<uint32_t>("maxSize") : 0) { }
00048             bool  isSet_, usesCharge_, usesSize_;
00049             float maxChi2_, maxCharge_;
00050             size_t maxSize_;
00051         };
00052         static const unsigned int NumberOfParamBlocks = 6;
00053 
00054         edm::InputTag trajectories_;
00055         bool doStrip_, doPixel_;
00056         edm::InputTag stripClusters_, pixelClusters_;
00057         bool mergeOld_;
00058         edm::InputTag oldRemovalInfo_;
00059 
00060         ParamBlock pblocks_[NumberOfParamBlocks];
00061         void readPSet(const edm::ParameterSet& iConfig, const std::string &name, 
00062                 int id1=-1, int id2=-1, int id3=-1, int id4=-1, int id5=-1, int id6=-1) ;
00063 
00064         std::vector<uint8_t> pixels, strips;                // avoid unneed alloc/dealloc of this
00065         edm::ProductID pixelSourceProdID, stripSourceProdID; // ProdIDs refs must point to (for consistency tests)
00066 
00067         inline void process(const TrackingRecHit *hit, float chi2);
00068         inline void process(const SiStripRecHit2D *hit2D, uint32_t subdet);
00069         inline void process(const SiStripRecHit1D *hit1D, 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 
00080 
00081 using namespace std;
00082 using namespace edm;
00083 using namespace reco;
00084 
00085 void
00086 TrackClusterRemover::readPSet(const edm::ParameterSet& iConfig, const std::string &name, 
00087     int id1, int id2, int id3, int id4, int id5, int id6) 
00088 {
00089     if (iConfig.exists(name)) {
00090         ParamBlock pblock(iConfig.getParameter<ParameterSet>(name));
00091         if (id1 == -1) {
00092             fill(pblocks_, pblocks_+NumberOfParamBlocks, pblock);
00093         } else {
00094             pblocks_[id1] = pblock;
00095             if (id2 != -1) pblocks_[id2] = pblock;
00096             if (id3 != -1) pblocks_[id3] = pblock;
00097             if (id4 != -1) pblocks_[id4] = pblock;
00098             if (id5 != -1) pblocks_[id5] = pblock;
00099             if (id6 != -1) pblocks_[id6] = pblock;
00100         }
00101     }
00102 }
00103 
00104 TrackClusterRemover::TrackClusterRemover(const ParameterSet& iConfig):
00105     trajectories_(iConfig.getParameter<InputTag>("trajectories")),
00106     doStrip_(iConfig.existsAs<bool>("doStrip") ? iConfig.getParameter<bool>("doStrip") : true),
00107     doPixel_(iConfig.existsAs<bool>("doPixel") ? iConfig.getParameter<bool>("doPixel") : true),
00108     stripClusters_(doStrip_ ? iConfig.getParameter<InputTag>("stripClusters") : InputTag("NONE")),
00109     pixelClusters_(doPixel_ ? iConfig.getParameter<InputTag>("pixelClusters") : InputTag("NONE")),
00110     mergeOld_(iConfig.exists("oldClusterRemovalInfo")),
00111     oldRemovalInfo_(mergeOld_ ? iConfig.getParameter<InputTag>("oldClusterRemovalInfo") : InputTag("NONE"))
00112 {
00113     if (doPixel_) produces< edmNew::DetSetVector<SiPixelCluster> >();
00114     if (doStrip_) produces< edmNew::DetSetVector<SiStripCluster> >();
00115     produces< ClusterRemovalInfo >();
00116 
00117     fill(pblocks_, pblocks_+NumberOfParamBlocks, ParamBlock());
00118     readPSet(iConfig, "Common",-1);
00119     if (doPixel_) {
00120         readPSet(iConfig, "Pixel" ,0,1);
00121         readPSet(iConfig, "PXB" ,0);
00122         readPSet(iConfig, "PXE" ,1);
00123     }
00124     if (doStrip_) {
00125         readPSet(iConfig, "Strip" ,2,3,4,5);
00126         readPSet(iConfig, "StripInner" ,2,3);
00127         readPSet(iConfig, "StripOuter" ,4,5);
00128         readPSet(iConfig, "TIB" ,2);
00129         readPSet(iConfig, "TID" ,3);
00130         readPSet(iConfig, "TOB" ,4);
00131         readPSet(iConfig, "TEC" ,5);
00132     }
00133 
00134     bool usingCharge = false;
00135     for (size_t i = 0; i < NumberOfParamBlocks; ++i) {
00136         if (!pblocks_[i].isSet_) throw cms::Exception("Configuration Error") << "TrackClusterRemover: Missing configuration for detector with subDetID = " << (i+1);
00137         if (pblocks_[i].usesCharge_ && !usingCharge) {
00138             throw cms::Exception("Configuration Error") << "TrackClusterRemover: Configuration for subDetID = " << (i+1) << " uses cluster charge, which is not enabled.";
00139         }
00140     }
00141 }
00142 
00143 
00144 TrackClusterRemover::~TrackClusterRemover()
00145 {
00146 }
00147 
00148 void TrackClusterRemover::mergeOld(ClusterRemovalInfo::Indices &refs,
00149                                             const ClusterRemovalInfo::Indices &oldRefs) 
00150 {
00151         for (size_t i = 0, n = refs.size(); i < n; ++i) {
00152             refs[i] = oldRefs[refs[i]];
00153        }
00154 }
00155 
00156  
00157 template<typename T> 
00158 auto_ptr<edmNew::DetSetVector<T> >
00159 TrackClusterRemover::cleanup(const edmNew::DetSetVector<T> &oldClusters, const std::vector<uint8_t> &isGood, 
00160                                 reco::ClusterRemovalInfo::Indices &refs, const reco::ClusterRemovalInfo::Indices *oldRefs) {  
00161     typedef typename edmNew::DetSetVector<T>             DSV;
00162     typedef typename edmNew::DetSetVector<T>::FastFiller DSF;
00163     typedef typename edmNew::DetSet<T>                   DS;
00164     auto_ptr<DSV> output(new DSV());
00165     output->reserve(oldClusters.size(), oldClusters.dataSize());
00166 
00167     // cluster removal loop
00168     const T * firstOffset = & oldClusters.data().front();
00169     for (typename DSV::const_iterator itdet = oldClusters.begin(), enddet = oldClusters.end(); itdet != enddet; ++itdet) {
00170         DS oldDS = *itdet;
00171 
00172         if (oldDS.empty()) continue; // skip empty detsets
00173 
00174         uint32_t id = oldDS.detId();
00175         DSF outds(*output, id);
00176 
00177         for (typename DS::const_iterator it = oldDS.begin(), ed = oldDS.end(); it != ed; ++it) {
00178             uint32_t index = ((&*it) - firstOffset);
00179             if (isGood[index]) { 
00180                 outds.push_back(*it);
00181                 refs.push_back(index); 
00182                 //std::cout << "TrackClusterRemover::cleanup " << typeid(T).name() << " reference " << index << " to " << (refs.size() - 1) << std::endl;
00183             }
00184         }
00185         if (outds.empty()) outds.abort(); // not write in an empty DSV
00186     }
00187     if (oldRefs != 0) mergeOld(refs, *oldRefs);
00188     return output;
00189 }
00190 
00191 void TrackClusterRemover::process(const SiStripRecHit2D *hit, uint32_t subdet) {
00192     SiStripRecHit2D::ClusterRef cluster = hit->cluster();
00193 
00194     if (cluster.id() != stripSourceProdID) throw cms::Exception("Inconsistent Data") << 
00195             "TrackClusterRemover: strip cluster ref from Product ID = " << cluster.id() << 
00196             " does not match with source cluster collection (ID = " << stripSourceProdID << ")\n.";
00197 
00198     assert(cluster.id() == stripSourceProdID);
00199     if (pblocks_[subdet-1].usesSize_ && (cluster->amplitudes().size() > pblocks_[subdet-1].maxSize_)) return;
00200 
00201 //DBG// cout << "Individual HIT " << cluster.key().first << ", INDEX = " << cluster.key().second << endl;
00202     strips[cluster.key()] = false;
00203 }
00204 
00205 void TrackClusterRemover::process(const SiStripRecHit1D *hit, uint32_t subdet) {
00206     SiStripRecHit1D::ClusterRef cluster = hit->cluster();
00207 
00208     if (cluster.id() != stripSourceProdID) throw cms::Exception("Inconsistent Data") << 
00209             "TrackClusterRemover: strip cluster ref from Product ID = " << cluster.id() << 
00210             " does not match with source cluster collection (ID = " << stripSourceProdID << ")\n.";
00211 
00212     assert(cluster.id() == stripSourceProdID);
00213     if (pblocks_[subdet-1].usesSize_ && (cluster->amplitudes().size() > pblocks_[subdet-1].maxSize_)) return;
00214 
00215 //DBG// cout << "Individual HIT " << cluster.key().first << ", INDEX = " << cluster.key().second << endl;
00216     strips[cluster.key()] = false;
00217 }
00218 
00219 
00220 
00221 void TrackClusterRemover::process(const TrackingRecHit *hit, float chi2) {
00222     DetId detid = hit->geographicalId(); 
00223     uint32_t subdet = detid.subdetId();
00224 
00225     assert ((subdet > 0) && (subdet <= NumberOfParamBlocks));
00226 
00227     // chi2 cut
00228     if (chi2 > pblocks_[subdet-1].maxChi2_) return;
00229 
00230     if ((subdet == PixelSubdetector::PixelBarrel) || (subdet == PixelSubdetector::PixelEndcap)) {
00231         if (!doPixel_) return;
00232         // this is a pixel, and i *know* it is
00233         const SiPixelRecHit *pixelHit = static_cast<const SiPixelRecHit *>(hit);
00234 
00235         SiPixelRecHit::ClusterRef cluster = pixelHit->cluster();
00236 
00237         if (cluster.id() != pixelSourceProdID) throw cms::Exception("Inconsistent Data") << 
00238                 "TrackClusterRemover: pixel cluster ref from Product ID = " << cluster.id() << 
00239                 " does not match with source cluster collection (ID = " << pixelSourceProdID << ")\n.";
00240 
00241         assert(cluster.id() == pixelSourceProdID);
00242 //DBG// cout << "HIT NEW PIXEL DETID = " << detid.rawId() << ", Cluster [ " << cluster.key().first << " / " <<  cluster.key().second << " ] " << endl;
00243 
00244         // if requested, cut on cluster size
00245         if (pblocks_[subdet-1].usesSize_ && (cluster->pixels().size() > pblocks_[subdet-1].maxSize_)) return;
00246 
00247         // mark as used
00248         pixels[cluster.key()] = false;
00249     } else { // aka Strip
00250         if (!doStrip_) return;
00251         const type_info &hitType = typeid(*hit);
00252         if (hitType == typeid(SiStripRecHit2D)) {
00253             const SiStripRecHit2D *stripHit = static_cast<const SiStripRecHit2D *>(hit);
00254 //DBG//     cout << "Plain RecHit 2D: " << endl;
00255             process(stripHit,subdet);}
00256         else if (hitType == typeid(SiStripRecHit1D)) {
00257           const SiStripRecHit1D *hit1D = static_cast<const SiStripRecHit1D *>(hit);
00258           process(hit1D,subdet);
00259         } else if (hitType == typeid(SiStripMatchedRecHit2D)) {
00260             const SiStripMatchedRecHit2D *matchHit = static_cast<const SiStripMatchedRecHit2D *>(hit);
00261 //DBG//     cout << "Matched RecHit 2D: " << endl;
00262             process(matchHit->monoHit(),subdet);
00263             process(matchHit->stereoHit(),subdet);
00264         } else if (hitType == typeid(ProjectedSiStripRecHit2D)) {
00265             const ProjectedSiStripRecHit2D *projHit = static_cast<const ProjectedSiStripRecHit2D *>(hit);
00266 //DBG//     cout << "Projected RecHit 2D: " << endl;
00267             process(&projHit->originalHit(),subdet);
00268         } else throw cms::Exception("NOT IMPLEMENTED") << "Don't know how to handle " << hitType.name() << " on detid " << detid.rawId() << "\n";
00269     }
00270 }
00271 
00272 /*   Schematic picture of n-th step Iterative removal
00273  *   (that os removing clusters after n-th step tracking)
00274  *   clusters:    [ C1 ] -> [ C2 ] -> ... -> [ Cn ] -> [ Cn + 1 ]
00275  *                   ^                          ^           ^--- OUTPUT "new" ID
00276  *                   |-- before any removal     |----- Source clusters                                                                             
00277  *                   |-- OUTPUT "old" ID        |----- Hits in Traj. point here                       
00278  *                   |                          \----- Old ClusterRemovalInfo "new" ID                 
00279  *                   \-- Old ClusterRemovalInfo "old" ID                                                  
00280  */
00281 
00282 
00283 void
00284 TrackClusterRemover::produce(Event& iEvent, const EventSetup& iSetup)
00285 {
00286     ProductID pixelOldProdID, stripOldProdID;
00287 
00288     Handle<edmNew::DetSetVector<SiPixelCluster> > pixelClusters;
00289     if (doPixel_) {
00290         iEvent.getByLabel(pixelClusters_, pixelClusters);
00291         pixelSourceProdID = pixelClusters.id();
00292     }
00293 //DBG// std::cout << "TrackClusterRemover: Read pixel " << pixelClusters_.encode() << " = ID " << pixelSourceProdID << std::endl;
00294 
00295     Handle<edmNew::DetSetVector<SiStripCluster> > stripClusters;
00296     if (doStrip_) {
00297         iEvent.getByLabel(stripClusters_, stripClusters);
00298         stripSourceProdID = stripClusters.id();
00299     }
00300 //DBG// std::cout << "TrackClusterRemover: Read strip " << stripClusters_.encode() << " = ID " << stripSourceProdID << std::endl;
00301 
00302     auto_ptr<ClusterRemovalInfo> cri;
00303     if (doStrip_ && doPixel_) cri.reset(new ClusterRemovalInfo(pixelClusters, stripClusters));
00304     else if (doStrip_) cri.reset(new ClusterRemovalInfo(stripClusters));
00305     else if (doPixel_) cri.reset(new ClusterRemovalInfo(pixelClusters));
00306 
00307     Handle<ClusterRemovalInfo> oldRemovalInfo;
00308     if (mergeOld_) { 
00309         iEvent.getByLabel(oldRemovalInfo_, oldRemovalInfo); 
00310         // Check ProductIDs
00311         if ( (oldRemovalInfo->stripNewRefProd().id() == stripClusters.id()) &&
00312              (oldRemovalInfo->pixelNewRefProd().id() == pixelClusters.id()) ) {
00313 
00314             cri->getOldClustersFrom(*oldRemovalInfo);
00315 
00316             pixelOldProdID = oldRemovalInfo->pixelRefProd().id();
00317             stripOldProdID = oldRemovalInfo->stripRefProd().id();
00318 
00319         } else {
00320             throw cms::Exception("Inconsistent Data") << "TrackClusterRemover: " <<
00321                 "Input collection product IDs are [pixel: " << pixelClusters.id() << ", strip: " << stripClusters.id() << "] \n" <<
00322                 "\t but the *old* ClusterRemovalInfo " << oldRemovalInfo_.encode() << " refers as 'new product ids' to " <<
00323                     "[pixel: " << oldRemovalInfo->pixelNewRefProd().id() << ", strip: " << oldRemovalInfo->stripNewRefProd().id() << "]\n" << 
00324                 "NOTA BENE: when running TrackClusterRemover with an old ClusterRemovalInfo the hits in the trajectory MUST be already re-keyed.\n";
00325         }
00326     } else { // then Old == Source
00327         pixelOldProdID = pixelSourceProdID;
00328         stripOldProdID = stripSourceProdID;
00329     }
00330 
00331     //Handle<TrajTrackAssociationCollection> trajectories; 
00332     Handle<vector<Trajectory> > trajectories; 
00333     iEvent.getByLabel(trajectories_, trajectories);
00334 
00335     if (doStrip_) {
00336         strips.resize(stripClusters->dataSize()); fill(strips.begin(), strips.end(), true);
00337     }
00338     if (doPixel_) {
00339         pixels.resize(pixelClusters->dataSize()); fill(pixels.begin(), pixels.end(), true);
00340     }
00341 
00342     //for (TrajTrackAssociationCollection::const_iterator it = trajectories->begin(), ed = trajectories->end(); it != ed; ++it)  {
00343     //    const Trajectory &tj = * it->key;
00344     for (vector<Trajectory>::const_iterator it = trajectories->begin(), ed = trajectories->end(); it != ed; ++it)  {
00345         const Trajectory &tj = * it;
00346         const vector<TrajectoryMeasurement> &tms = tj.measurements();
00347         vector<TrajectoryMeasurement>::const_iterator itm, endtm;
00348         for (itm = tms.begin(), endtm = tms.end(); itm != endtm; ++itm) {
00349             const TrackingRecHit *hit = itm->recHit()->hit();
00350             if (!hit->isValid()) continue; 
00351             process( hit, itm->estimate() );
00352         }
00353     }
00354 
00355     if (doPixel_) {
00356         auto_ptr<edmNew::DetSetVector<SiPixelCluster> > newPixelClusters = cleanup(*pixelClusters, pixels, 
00357                     cri->pixelIndices(), mergeOld_ ? &oldRemovalInfo->pixelIndices() : 0);
00358         OrphanHandle<edmNew::DetSetVector<SiPixelCluster> > newPixels = iEvent.put(newPixelClusters); 
00359 //DBG// std::cout << "TrackClusterRemover: Wrote pixel " << newPixels.id() << " from " << pixelSourceProdID << std::endl;
00360         cri->setNewPixelClusters(newPixels);
00361     }
00362     if (doStrip_) {
00363         auto_ptr<edmNew::DetSetVector<SiStripCluster> > newStripClusters = cleanup(*stripClusters, strips, 
00364                     cri->stripIndices(), mergeOld_ ? &oldRemovalInfo->stripIndices() : 0);
00365         OrphanHandle<edmNew::DetSetVector<SiStripCluster> > newStrips = iEvent.put(newStripClusters); 
00366 //DBG// std::cout << "TrackClusterRemover: Wrote strip " << newStrips.id() << " from " << stripSourceProdID << std::endl;
00367         cri->setNewStripClusters(newStrips);
00368     }
00369 
00370 
00371     iEvent.put(cri);
00372 
00373     pixels.clear(); strips.clear(); 
00374 }
00375 
00376 #include "FWCore/PluginManager/interface/ModuleDef.h"
00377 #include "FWCore/Framework/interface/MakerMacros.h"
00378 DEFINE_FWK_MODULE(TrackClusterRemover);