Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2021-03-14 23:24:27

0001 #include "FWCore/Framework/interface/Frameworkfwd.h"
0002 #include "FWCore/Framework/interface/stream/EDProducer.h"
0003 #include "FWCore/Framework/interface/Event.h"
0004 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0005 #include "FWCore/Utilities/interface/InputTag.h"
0006 
0007 #include "DataFormats/SiStripCluster/interface/SiStripCluster.h"
0008 #include "DataFormats/SiPixelCluster/interface/SiPixelCluster.h"
0009 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2D.h"
0010 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit1D.h"
0011 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2D.h"
0012 #include "DataFormats/TrackerRecHit2D/interface/ProjectedSiStripRecHit2D.h"
0013 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
0014 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHitCollection.h"
0015 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2DCollection.h"
0016 #include "DataFormats/Common/interface/Handle.h"
0017 #include "DataFormats/SiStripDetId/interface/SiStripDetId.h"
0018 #include "DataFormats/Common/interface/DetSetVector.h"
0019 #include "DataFormats/Common/interface/ValueMap.h"
0020 #include "DataFormats/Common/interface/DetSetVectorNew.h"
0021 #include "DataFormats/Provenance/interface/ProductID.h"
0022 #include "DataFormats/Common/interface/ContainerMask.h"
0023 
0024 #include "DataFormats/DetId/interface/DetId.h"
0025 
0026 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0027 #include "DataFormats/TrackReco/interface/Track.h"
0028 #include "DataFormats/TrackerRecHit2D/interface/ClusterRemovalInfo.h"
0029 
0030 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0031 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0032 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
0033 #include "Geometry/CommonDetUnit/interface/GeomDetType.h"
0034 
0035 #include "RecoTracker/TransientTrackingRecHit/interface/Traj2TrackHits.h"
0036 
0037 //
0038 // class decleration
0039 //
0040 
0041 class HITrackClusterRemover : public edm::stream::EDProducer<> {
0042 public:
0043   HITrackClusterRemover(const edm::ParameterSet &iConfig);
0044   ~HITrackClusterRemover() override;
0045   void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override;
0046 
0047 private:
0048   edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> const tTrackerGeom_;
0049   struct ParamBlock {
0050     ParamBlock() : isSet_(false), usesCharge_(false) {}
0051     ParamBlock(const edm::ParameterSet &iConfig)
0052         : isSet_(true),
0053           usesCharge_(iConfig.exists("maxCharge")),
0054           usesSize_(iConfig.exists("maxSize")),
0055           cutOnPixelCharge_(iConfig.exists("minGoodPixelCharge")),
0056           cutOnStripCharge_(iConfig.exists("minGoodStripCharge")),
0057           maxChi2_(iConfig.getParameter<double>("maxChi2")),
0058           maxCharge_(usesCharge_ ? iConfig.getParameter<double>("maxCharge") : 0),
0059           minGoodPixelCharge_(cutOnPixelCharge_ ? iConfig.getParameter<double>("minGoodPixelCharge") : 0),
0060           minGoodStripCharge_(cutOnStripCharge_ ? iConfig.getParameter<double>("minGoodStripCharge") : 0),
0061           maxSize_(usesSize_ ? iConfig.getParameter<uint32_t>("maxSize") : 0) {}
0062     bool isSet_, usesCharge_, usesSize_, cutOnPixelCharge_, cutOnStripCharge_;
0063     float maxChi2_, maxCharge_, minGoodPixelCharge_, minGoodStripCharge_;
0064     size_t maxSize_;
0065   };
0066   static const unsigned int NumberOfParamBlocks = 6;
0067 
0068   bool doTracks_;
0069   bool doStrip_, doPixel_;
0070   bool mergeOld_;
0071   typedef edm::ContainerMask<edmNew::DetSetVector<SiPixelCluster> > PixelMaskContainer;
0072   typedef edm::ContainerMask<edmNew::DetSetVector<SiStripCluster> > StripMaskContainer;
0073   edm::EDGetTokenT<edmNew::DetSetVector<SiPixelCluster> > pixelClusters_;
0074   edm::EDGetTokenT<edmNew::DetSetVector<SiStripCluster> > stripClusters_;
0075   edm::EDGetTokenT<reco::TrackCollection> tracks_;
0076   edm::EDGetTokenT<reco::ClusterRemovalInfo> oldRemovalInfo_;
0077   edm::EDGetTokenT<PixelMaskContainer> oldPxlMaskToken_;
0078   edm::EDGetTokenT<StripMaskContainer> oldStrMaskToken_;
0079   std::vector<edm::EDGetTokenT<edm::ValueMap<int> > > overrideTrkQuals_;
0080   edm::EDGetTokenT<SiStripRecHit2DCollection> rphiRecHitToken_, stereoRecHitToken_;
0081   //    edm::EDGetTokenT<SiPixelRecHitCollection> pixelRecHitsToken_;
0082 
0083   ParamBlock pblocks_[NumberOfParamBlocks];
0084   void readPSet(const edm::ParameterSet &iConfig,
0085                 const std::string &name,
0086                 int id1 = -1,
0087                 int id2 = -1,
0088                 int id3 = -1,
0089                 int id4 = -1,
0090                 int id5 = -1,
0091                 int id6 = -1);
0092 
0093   std::vector<uint8_t> pixels, strips;                  // avoid unneed alloc/dealloc of this
0094   edm::ProductID pixelSourceProdID, stripSourceProdID;  // ProdIDs refs must point to (for consistency tests)
0095 
0096   inline void process(const TrackingRecHit *hit, unsigned char chi2, const TrackerGeometry *tg);
0097   inline void process(const OmniClusterRef &cluRef, SiStripDetId &detid, bool fromTrack);
0098 
0099   template <typename T>
0100   std::unique_ptr<edmNew::DetSetVector<T> > cleanup(const edmNew::DetSetVector<T> &oldClusters,
0101                                                     const std::vector<uint8_t> &isGood,
0102                                                     reco::ClusterRemovalInfo::Indices &refs,
0103                                                     const reco::ClusterRemovalInfo::Indices *oldRefs);
0104 
0105   // Carries in full removal info about a given det from oldRefs
0106   void mergeOld(reco::ClusterRemovalInfo::Indices &refs, const reco::ClusterRemovalInfo::Indices &oldRefs);
0107 
0108   bool clusterWasteSolution_, doStripChargeCheck_, doPixelChargeCheck_;
0109   std::string stripRecHits_, pixelRecHits_;
0110   bool filterTracks_;
0111   int minNumberOfLayersWithMeasBeforeFiltering_;
0112   reco::TrackBase::TrackQuality trackQuality_;
0113   std::vector<bool> collectedStrips_;
0114   std::vector<bool> collectedPixels_;
0115 
0116   float sensorThickness(const SiStripDetId &detid) const;
0117 };
0118 
0119 using namespace std;
0120 using namespace edm;
0121 using namespace reco;
0122 
0123 void HITrackClusterRemover::readPSet(
0124     const edm::ParameterSet &iConfig, const std::string &name, int id1, int id2, int id3, int id4, int id5, int id6) {
0125   if (iConfig.exists(name)) {
0126     ParamBlock pblock(iConfig.getParameter<ParameterSet>(name));
0127     if (id1 == -1) {
0128       fill(pblocks_, pblocks_ + NumberOfParamBlocks, pblock);
0129     } else {
0130       pblocks_[id1] = pblock;
0131       if (id2 != -1)
0132         pblocks_[id2] = pblock;
0133       if (id3 != -1)
0134         pblocks_[id3] = pblock;
0135       if (id4 != -1)
0136         pblocks_[id4] = pblock;
0137       if (id5 != -1)
0138         pblocks_[id5] = pblock;
0139       if (id6 != -1)
0140         pblocks_[id6] = pblock;
0141     }
0142   }
0143 }
0144 
0145 HITrackClusterRemover::HITrackClusterRemover(const ParameterSet &iConfig)
0146     : tTrackerGeom_(esConsumes<TrackerGeometry, TrackerDigiGeometryRecord>()),
0147       doTracks_(iConfig.exists("trajectories")),
0148       doStrip_(iConfig.existsAs<bool>("doStrip") ? iConfig.getParameter<bool>("doStrip") : true),
0149       doPixel_(iConfig.existsAs<bool>("doPixel") ? iConfig.getParameter<bool>("doPixel") : true),
0150       mergeOld_(iConfig.exists("oldClusterRemovalInfo")),
0151       clusterWasteSolution_(true),
0152       doStripChargeCheck_(
0153           iConfig.existsAs<bool>("doStripChargeCheck") ? iConfig.getParameter<bool>("doStripChargeCheck") : false),
0154       doPixelChargeCheck_(
0155           iConfig.existsAs<bool>("doPixelChargeCheck") ? iConfig.getParameter<bool>("doPixelChargeCheck") : false),
0156       stripRecHits_(doStripChargeCheck_ ? iConfig.getParameter<std::string>("stripRecHits")
0157                                         : std::string("siStripMatchedRecHits")),
0158       pixelRecHits_(doPixelChargeCheck_ ? iConfig.getParameter<std::string>("pixelRecHits")
0159                                         : std::string("siPixelRecHits")) {
0160   mergeOld_ = mergeOld_ && !iConfig.getParameter<InputTag>("oldClusterRemovalInfo").label().empty();
0161   if (iConfig.exists("overrideTrkQuals"))
0162     overrideTrkQuals_.push_back(consumes<edm::ValueMap<int> >(iConfig.getParameter<InputTag>("overrideTrkQuals")));
0163   if (iConfig.exists("clusterLessSolution"))
0164     clusterWasteSolution_ = !iConfig.getParameter<bool>("clusterLessSolution");
0165   if ((doPixelChargeCheck_ && !doPixel_) || (doStripChargeCheck_ && !doStrip_))
0166     throw cms::Exception("Configuration Error")
0167         << "HITrackClusterRemover: Charge check asked without cluster collection ";
0168   if (doPixelChargeCheck_)
0169     throw cms::Exception("Configuration Error")
0170         << "HITrackClusterRemover: Pixel cluster charge check not yet implemented";
0171 
0172   if (doPixel_ && clusterWasteSolution_)
0173     produces<edmNew::DetSetVector<SiPixelCluster> >();
0174   if (doStrip_ && clusterWasteSolution_)
0175     produces<edmNew::DetSetVector<SiStripCluster> >();
0176   if (clusterWasteSolution_)
0177     produces<ClusterRemovalInfo>();
0178 
0179   assert(!clusterWasteSolution_);
0180 
0181   fill(pblocks_, pblocks_ + NumberOfParamBlocks, ParamBlock());
0182   readPSet(iConfig, "Common", -1);
0183   if (doPixel_) {
0184     readPSet(iConfig, "Pixel", 0, 1);
0185     readPSet(iConfig, "PXB", 0);
0186     readPSet(iConfig, "PXE", 1);
0187   }
0188   if (doStrip_) {
0189     readPSet(iConfig, "Strip", 2, 3, 4, 5);
0190     readPSet(iConfig, "StripInner", 2, 3);
0191     readPSet(iConfig, "StripOuter", 4, 5);
0192     readPSet(iConfig, "TIB", 2);
0193     readPSet(iConfig, "TID", 3);
0194     readPSet(iConfig, "TOB", 4);
0195     readPSet(iConfig, "TEC", 5);
0196   }
0197 
0198   bool usingCharge = false;
0199   for (size_t i = 0; i < NumberOfParamBlocks; ++i) {
0200     if (!pblocks_[i].isSet_)
0201       throw cms::Exception("Configuration Error")
0202           << "HITrackClusterRemover: Missing configuration for detector with subDetID = " << (i + 1);
0203     if (pblocks_[i].usesCharge_ && !usingCharge) {
0204       throw cms::Exception("Configuration Error") << "HITrackClusterRemover: Configuration for subDetID = " << (i + 1)
0205                                                   << " uses cluster charge, which is not enabled.";
0206     }
0207   }
0208 
0209   if (!clusterWasteSolution_) {
0210     produces<edm::ContainerMask<edmNew::DetSetVector<SiPixelCluster> > >();
0211     produces<edm::ContainerMask<edmNew::DetSetVector<SiStripCluster> > >();
0212   }
0213   trackQuality_ = reco::TrackBase::undefQuality;
0214   filterTracks_ = false;
0215   if (iConfig.exists("TrackQuality")) {
0216     filterTracks_ = true;
0217     trackQuality_ = reco::TrackBase::qualityByName(iConfig.getParameter<std::string>("TrackQuality"));
0218     minNumberOfLayersWithMeasBeforeFiltering_ =
0219         iConfig.existsAs<int>("minNumberOfLayersWithMeasBeforeFiltering")
0220             ? iConfig.getParameter<int>("minNumberOfLayersWithMeasBeforeFiltering")
0221             : 0;
0222   }
0223 
0224   if (doTracks_)
0225     tracks_ = consumes<reco::TrackCollection>(iConfig.getParameter<InputTag>("trajectories"));
0226   if (doPixel_)
0227     pixelClusters_ = consumes<edmNew::DetSetVector<SiPixelCluster> >(iConfig.getParameter<InputTag>("pixelClusters"));
0228   if (doStrip_)
0229     stripClusters_ = consumes<edmNew::DetSetVector<SiStripCluster> >(iConfig.getParameter<InputTag>("stripClusters"));
0230   if (mergeOld_) {
0231     oldRemovalInfo_ = consumes<ClusterRemovalInfo>(iConfig.getParameter<InputTag>("oldClusterRemovalInfo"));
0232     oldPxlMaskToken_ = consumes<PixelMaskContainer>(iConfig.getParameter<InputTag>("oldClusterRemovalInfo"));
0233     oldStrMaskToken_ = consumes<StripMaskContainer>(iConfig.getParameter<InputTag>("oldClusterRemovalInfo"));
0234   }
0235 
0236   if (doStripChargeCheck_) {
0237     rphiRecHitToken_ = consumes<SiStripRecHit2DCollection>(InputTag(stripRecHits_, "rphiRecHit"));
0238     stereoRecHitToken_ = consumes<SiStripRecHit2DCollection>(InputTag(stripRecHits_, "stereoRecHit"));
0239   }
0240   //    if(doPixelChargeCheck_) pixelRecHitsToken_ = consumes<SiPixelRecHitCollection>(InputTag(pixelRecHits_));
0241 }
0242 
0243 HITrackClusterRemover::~HITrackClusterRemover() {}
0244 
0245 void HITrackClusterRemover::mergeOld(ClusterRemovalInfo::Indices &refs, const ClusterRemovalInfo::Indices &oldRefs) {
0246   for (size_t i = 0, n = refs.size(); i < n; ++i) {
0247     refs[i] = oldRefs[refs[i]];
0248   }
0249 }
0250 
0251 template <typename T>
0252 std::unique_ptr<edmNew::DetSetVector<T> > HITrackClusterRemover::cleanup(
0253     const edmNew::DetSetVector<T> &oldClusters,
0254     const std::vector<uint8_t> &isGood,
0255     reco::ClusterRemovalInfo::Indices &refs,
0256     const reco::ClusterRemovalInfo::Indices *oldRefs) {
0257   typedef typename edmNew::DetSetVector<T> DSV;
0258   typedef typename edmNew::DetSetVector<T>::FastFiller DSF;
0259   typedef typename edmNew::DetSet<T> DS;
0260   auto output = std::make_unique<DSV>();
0261   output->reserve(oldClusters.size(), oldClusters.dataSize());
0262 
0263   unsigned int countOld = 0;
0264   unsigned int countNew = 0;
0265 
0266   // cluster removal loop
0267   const T *firstOffset = &oldClusters.data().front();
0268   for (typename DSV::const_iterator itdet = oldClusters.begin(), enddet = oldClusters.end(); itdet != enddet; ++itdet) {
0269     DS oldDS = *itdet;
0270 
0271     if (oldDS.empty())
0272       continue;  // skip empty detsets
0273 
0274     uint32_t id = oldDS.detId();
0275     DSF outds(*output, id);
0276 
0277     for (typename DS::const_iterator it = oldDS.begin(), ed = oldDS.end(); it != ed; ++it) {
0278       uint32_t index = ((&*it) - firstOffset);
0279       countOld++;
0280       if (isGood[index]) {
0281         outds.push_back(*it);
0282         countNew++;
0283         refs.push_back(index);
0284         //std::cout << "HITrackClusterRemover::cleanup " << typeid(T).name() << " reference " << index << " to " << (refs.size() - 1) << std::endl;
0285       }
0286     }
0287     if (outds.empty())
0288       outds.abort();  // not write in an empty DSV
0289   }
0290 
0291   if (oldRefs != nullptr)
0292     mergeOld(refs, *oldRefs);
0293   return output;
0294 }
0295 
0296 float HITrackClusterRemover::sensorThickness(const SiStripDetId &detid) const {
0297   if (detid.subdetId() >= SiStripDetId::TIB) {
0298     if (detid.subdetId() == SiStripDetId::TOB)
0299       return 0.047;
0300     if (detid.moduleGeometry() == SiStripModuleGeometry::W5 || detid.moduleGeometry() == SiStripModuleGeometry::W6 ||
0301         detid.moduleGeometry() == SiStripModuleGeometry::W7)
0302       return 0.047;
0303     return 0.029;  // so it is TEC ring 1-4 or TIB or TOB;
0304   } else if (detid.subdetId() == PixelSubdetector::PixelBarrel)
0305     return 0.0285;
0306   else
0307     return 0.027;
0308 }
0309 
0310 void HITrackClusterRemover::process(OmniClusterRef const &ocluster, SiStripDetId &detid, bool fromTrack) {
0311   SiStripRecHit2D::ClusterRef cluster = ocluster.cluster_strip();
0312   if (cluster.id() != stripSourceProdID)
0313     throw cms::Exception("Inconsistent Data")
0314         << "HITrackClusterRemover: strip cluster ref from Product ID = " << cluster.id()
0315         << " does not match with source cluster collection (ID = " << stripSourceProdID << ")\n.";
0316 
0317   uint32_t subdet = detid.subdetId();
0318   assert(cluster.id() == stripSourceProdID);
0319   if (pblocks_[subdet - 1].usesSize_ && (cluster->amplitudes().size() > pblocks_[subdet - 1].maxSize_))
0320     return;
0321   if (!fromTrack) {
0322     if (pblocks_[subdet - 1].cutOnStripCharge_ &&
0323         (cluster->charge() > (pblocks_[subdet - 1].minGoodStripCharge_ * sensorThickness(detid))))
0324       return;
0325   }
0326 
0327   if (collectedStrips_.size() <= cluster.key())
0328     edm::LogError("BadCollectionSize") << collectedStrips_.size() << " is smaller than " << cluster.key();
0329 
0330   assert(collectedStrips_.size() > cluster.key());
0331   strips[cluster.key()] = false;
0332   if (!clusterWasteSolution_)
0333     collectedStrips_[cluster.key()] = true;
0334 }
0335 
0336 void HITrackClusterRemover::process(const TrackingRecHit *hit, unsigned char chi2, const TrackerGeometry *tg) {
0337   SiStripDetId detid = hit->geographicalId();
0338   uint32_t subdet = detid.subdetId();
0339 
0340   assert((subdet > 0) && (subdet <= NumberOfParamBlocks));
0341 
0342   // chi2 cut
0343   if (chi2 > Traj2TrackHits::toChi2x5(pblocks_[subdet - 1].maxChi2_))
0344     return;
0345 
0346   if (GeomDetEnumerators::isTrackerPixel(tg->geomDetSubDetector(subdet))) {
0347     if (!doPixel_)
0348       return;
0349     // this is a pixel, and i *know* it is
0350     const SiPixelRecHit *pixelHit = static_cast<const SiPixelRecHit *>(hit);
0351 
0352     SiPixelRecHit::ClusterRef cluster = pixelHit->cluster();
0353 
0354     if (cluster.id() != pixelSourceProdID)
0355       throw cms::Exception("Inconsistent Data")
0356           << "HITrackClusterRemover: pixel cluster ref from Product ID = " << cluster.id()
0357           << " does not match with source cluster collection (ID = " << pixelSourceProdID << ")\n.";
0358 
0359     assert(cluster.id() == pixelSourceProdID);
0360     //DBG// cout << "HIT NEW PIXEL DETID = " << detid.rawId() << ", Cluster [ " << cluster.key().first << " / " <<  cluster.key().second << " ] " << endl;
0361 
0362     // if requested, cut on cluster size
0363     if (pblocks_[subdet - 1].usesSize_ && (cluster->pixels().size() > pblocks_[subdet - 1].maxSize_))
0364       return;
0365 
0366     // mark as used
0367     pixels[cluster.key()] = false;
0368 
0369     //if(!clusterWasteSolution_) collectedPixel[detid.rawId()].insert(cluster);
0370     assert(collectedPixels_.size() > cluster.key());
0371     //assert(detid.rawId() == cluster->geographicalId()); //This condition fails
0372     if (!clusterWasteSolution_)
0373       collectedPixels_[cluster.key()] = true;
0374 
0375   } else {  // aka Strip
0376     if (!doStrip_)
0377       return;
0378     const type_info &hitType = typeid(*hit);
0379     if (hitType == typeid(SiStripRecHit2D)) {
0380       const SiStripRecHit2D *stripHit = static_cast<const SiStripRecHit2D *>(hit);
0381       //DBG//     cout << "Plain RecHit 2D: " << endl;
0382       process(stripHit->omniClusterRef(), detid, true);
0383     } else if (hitType == typeid(SiStripRecHit1D)) {
0384       const SiStripRecHit1D *hit1D = static_cast<const SiStripRecHit1D *>(hit);
0385       process(hit1D->omniClusterRef(), detid, true);
0386     } else if (hitType == typeid(SiStripMatchedRecHit2D)) {
0387       const SiStripMatchedRecHit2D *matchHit = static_cast<const SiStripMatchedRecHit2D *>(hit);
0388       //DBG//     cout << "Matched RecHit 2D: " << endl;
0389       process(matchHit->monoClusterRef(), detid, true);
0390       process(matchHit->stereoClusterRef(), detid, true);
0391     } else if (hitType == typeid(ProjectedSiStripRecHit2D)) {
0392       const ProjectedSiStripRecHit2D *projHit = static_cast<const ProjectedSiStripRecHit2D *>(hit);
0393       //DBG//     cout << "Projected RecHit 2D: " << endl;
0394       process(projHit->originalHit().omniClusterRef(), detid, true);
0395     } else
0396       throw cms::Exception("NOT IMPLEMENTED")
0397           << "Don't know how to handle " << hitType.name() << " on detid " << detid.rawId() << "\n";
0398   }
0399 }
0400 
0401 /*   Schematic picture of n-th step Iterative removal
0402  *   (that os removing clusters after n-th step tracking)
0403  *   clusters:    [ C1 ] -> [ C2 ] -> ... -> [ Cn ] -> [ Cn + 1 ]
0404  *                   ^                          ^           ^--- OUTPUT "new" ID
0405  *                   |-- before any removal     |----- Source clusters                                                                             
0406  *                   |-- OUTPUT "old" ID        |----- Hits in Traj. point here                       
0407  *                   |                          \----- Old ClusterRemovalInfo "new" ID                 
0408  *                   \-- Old ClusterRemovalInfo "old" ID                                                  
0409  */
0410 
0411 void HITrackClusterRemover::produce(Event &iEvent, const EventSetup &iSetup) {
0412   ProductID pixelOldProdID, stripOldProdID;
0413 
0414   const auto &tgh = &iSetup.getData(tTrackerGeom_);
0415 
0416   Handle<edmNew::DetSetVector<SiPixelCluster> > pixelClusters;
0417   if (doPixel_) {
0418     iEvent.getByToken(pixelClusters_, pixelClusters);
0419     pixelSourceProdID = pixelClusters.id();
0420   }
0421   //DBG// std::cout << "HITrackClusterRemover: Read pixel " << pixelClusters_.encode() << " = ID " << pixelSourceProdID << std::endl;
0422 
0423   Handle<edmNew::DetSetVector<SiStripCluster> > stripClusters;
0424   if (doStrip_) {
0425     iEvent.getByToken(stripClusters_, stripClusters);
0426     stripSourceProdID = stripClusters.id();
0427   }
0428   //DBG// std::cout << "HITrackClusterRemover: Read strip " << stripClusters_.encode() << " = ID " << stripSourceProdID << std::endl;
0429 
0430   std::unique_ptr<ClusterRemovalInfo> cri;
0431   if (clusterWasteSolution_) {
0432     if (doStrip_ && doPixel_)
0433       cri = std::make_unique<ClusterRemovalInfo>(pixelClusters, stripClusters);
0434     else if (doStrip_)
0435       cri = std::make_unique<ClusterRemovalInfo>(stripClusters);
0436     else if (doPixel_)
0437       cri = std::make_unique<ClusterRemovalInfo>(pixelClusters);
0438   }
0439 
0440   Handle<ClusterRemovalInfo> oldRemovalInfo;
0441   if (mergeOld_ && clusterWasteSolution_) {
0442     iEvent.getByToken(oldRemovalInfo_, oldRemovalInfo);
0443     // Check ProductIDs
0444     if ((oldRemovalInfo->stripNewRefProd().id() == stripClusters.id()) &&
0445         (oldRemovalInfo->pixelNewRefProd().id() == pixelClusters.id())) {
0446       cri->getOldClustersFrom(*oldRemovalInfo);
0447 
0448       pixelOldProdID = oldRemovalInfo->pixelRefProd().id();
0449       stripOldProdID = oldRemovalInfo->stripRefProd().id();
0450 
0451     } else {
0452       edm::EDConsumerBase::Labels labels;
0453       labelsForToken(oldRemovalInfo_, labels);
0454       throw cms::Exception("Inconsistent Data")
0455           << "HITrackClusterRemover: "
0456           << "Input collection product IDs are [pixel: " << pixelClusters.id() << ", strip: " << stripClusters.id()
0457           << "] \n"
0458           << "\t but the *old* ClusterRemovalInfo " << labels.productInstance << " refers as 'new product ids' to "
0459           << "[pixel: " << oldRemovalInfo->pixelNewRefProd().id()
0460           << ", strip: " << oldRemovalInfo->stripNewRefProd().id() << "]\n"
0461           << "NOTA BENE: when running HITrackClusterRemover with an old ClusterRemovalInfo the hits in the trajectory "
0462              "MUST be already re-keyed.\n";
0463     }
0464   } else {  // then Old == Source
0465     pixelOldProdID = pixelSourceProdID;
0466     stripOldProdID = stripSourceProdID;
0467   }
0468 
0469   if (doStrip_) {
0470     strips.resize(stripClusters->dataSize());
0471     fill(strips.begin(), strips.end(), true);
0472   }
0473   if (doPixel_) {
0474     pixels.resize(pixelClusters->dataSize());
0475     fill(pixels.begin(), pixels.end(), true);
0476   }
0477   if (mergeOld_) {
0478     edm::Handle<PixelMaskContainer> oldPxlMask;
0479     edm::Handle<StripMaskContainer> oldStrMask;
0480     iEvent.getByToken(oldPxlMaskToken_, oldPxlMask);
0481     iEvent.getByToken(oldStrMaskToken_, oldStrMask);
0482     LogDebug("HITrackClusterRemover") << "to merge in, " << oldStrMask->size() << " strp and " << oldPxlMask->size()
0483                                       << " pxl";
0484     oldStrMask->copyMaskTo(collectedStrips_);
0485     oldPxlMask->copyMaskTo(collectedPixels_);
0486     assert(stripClusters->dataSize() >= collectedStrips_.size());
0487     collectedStrips_.resize(stripClusters->dataSize(), false);
0488   } else {
0489     collectedStrips_.resize(stripClusters->dataSize(), false);
0490     collectedPixels_.resize(pixelClusters->dataSize(), false);
0491   }
0492 
0493   if (doTracks_) {
0494     Handle<reco::TrackCollection> tracks;
0495     iEvent.getByToken(tracks_, tracks);
0496 
0497     std::vector<Handle<edm::ValueMap<int> > > quals;
0498     if (!overrideTrkQuals_.empty()) {
0499       quals.resize(1);
0500       iEvent.getByToken(overrideTrkQuals_[0], quals[0]);
0501     }
0502     int it = 0;
0503     for (const auto &track : *tracks) {
0504       if (filterTracks_) {
0505         bool goodTk = true;
0506         if (!quals.empty()) {
0507           int qual = (*(quals[0])).get(it++);
0508           if (qual < 0) {
0509             goodTk = false;
0510           }
0511           //note that this does not work for some trackquals (goodIterative  or undefQuality)
0512           else
0513             goodTk = (qual & (1 << trackQuality_)) >> trackQuality_;
0514         } else
0515           goodTk = (track.quality(trackQuality_));
0516         if (!goodTk)
0517           continue;
0518         if (track.hitPattern().trackerLayersWithMeasurement() < minNumberOfLayersWithMeasBeforeFiltering_)
0519           continue;
0520       }
0521       auto const &chi2sX5 = track.extra()->chi2sX5();
0522       assert(chi2sX5.size() == track.recHitsSize());
0523       auto hb = track.recHitsBegin();
0524       for (unsigned int h = 0; h < track.recHitsSize(); h++) {
0525         auto hit = *(hb + h);
0526         if (!hit->isValid())
0527           continue;
0528         process(hit, chi2sX5[h], tgh);
0529       }
0530     }
0531   }
0532 
0533   if (doStripChargeCheck_) {
0534     edm::Handle<SiStripRecHit2DCollection> rechitsrphi;
0535     iEvent.getByToken(rphiRecHitToken_, rechitsrphi);
0536     const SiStripRecHit2DCollection::DataContainer *rphiRecHits = &(rechitsrphi).product()->data();
0537     for (SiStripRecHit2DCollection::DataContainer::const_iterator recHit = rphiRecHits->begin();
0538          recHit != rphiRecHits->end();
0539          recHit++) {
0540       SiStripDetId detid = recHit->geographicalId();
0541       process(recHit->omniClusterRef(), detid, false);
0542     }
0543     edm::Handle<SiStripRecHit2DCollection> rechitsstereo;
0544     iEvent.getByToken(stereoRecHitToken_, rechitsstereo);
0545     const SiStripRecHit2DCollection::DataContainer *stereoRecHits = &(rechitsstereo).product()->data();
0546     for (SiStripRecHit2DCollection::DataContainer::const_iterator recHit = stereoRecHits->begin();
0547          recHit != stereoRecHits->end();
0548          recHit++) {
0549       SiStripDetId detid = recHit->geographicalId();
0550       process(recHit->omniClusterRef(), detid, false);
0551     }
0552   }
0553   //    if(doPixelChargeCheck_) {
0554   //    edm::Handle<SiPixelRecHitCollection> pixelrechits;
0555   //    iEvent.getByToken(pixelRecHitsToken_,pixelrechits);
0556   //    }
0557 
0558   if (doPixel_ && clusterWasteSolution_) {
0559     OrphanHandle<edmNew::DetSetVector<SiPixelCluster> > newPixels = iEvent.put(
0560         cleanup(*pixelClusters, pixels, cri->pixelIndices(), mergeOld_ ? &oldRemovalInfo->pixelIndices() : nullptr));
0561     //DBG// std::cout << "HITrackClusterRemover: Wrote pixel " << newPixels.id() << " from " << pixelSourceProdID << std::endl;
0562     cri->setNewPixelClusters(newPixels);
0563   }
0564   if (doStrip_ && clusterWasteSolution_) {
0565     OrphanHandle<edmNew::DetSetVector<SiStripCluster> > newStrips = iEvent.put(
0566         cleanup(*stripClusters, strips, cri->stripIndices(), mergeOld_ ? &oldRemovalInfo->stripIndices() : nullptr));
0567     //DBG// std::cout << "HITrackClusterRemover: Wrote strip " << newStrips.id() << " from " << stripSourceProdID << std::endl;
0568     cri->setNewStripClusters(newStrips);
0569   }
0570 
0571   if (clusterWasteSolution_) {
0572     //      double fraction_pxl= cri->pixelIndices().size() / (double) pixels.size();
0573     //      double fraction_strp= cri->stripIndices().size() / (double) strips.size();
0574     //      edm::LogWarning("HITrackClusterRemover")<<" fraction: " << fraction_pxl <<" "<<fraction_strp;
0575     iEvent.put(std::move(cri));
0576   }
0577 
0578   pixels.clear();
0579   strips.clear();
0580 
0581   if (!clusterWasteSolution_) {
0582     //auto_ptr<edmNew::DetSetVector<SiPixelClusterRefNew> > removedPixelClsuterRefs(new edmNew::DetSetVector<SiPixelClusterRefNew>());
0583     //auto_ptr<edmNew::DetSetVector<SiStripRecHit1D::ClusterRef> > removedStripClsuterRefs(new edmNew::DetSetVector<SiStripRecHit1D::ClusterRef>());
0584 
0585     LogDebug("HITrackClusterRemover") << "total strip to skip: "
0586                                       << std::count(collectedStrips_.begin(), collectedStrips_.end(), true);
0587     // std::cout << "HITrackClusterRemover " <<"total strip to skip: "<<std::count(collectedStrips_.begin(),collectedStrips_.end(),true) <<std::endl;
0588     iEvent.put(std::make_unique<StripMaskContainer>(edm::RefProd<edmNew::DetSetVector<SiStripCluster> >(stripClusters),
0589                                                     collectedStrips_));
0590 
0591     LogDebug("HITrackClusterRemover") << "total pxl to skip: "
0592                                       << std::count(collectedPixels_.begin(), collectedPixels_.end(), true);
0593     iEvent.put(std::make_unique<PixelMaskContainer>(edm::RefProd<edmNew::DetSetVector<SiPixelCluster> >(pixelClusters),
0594                                                     collectedPixels_));
0595   }
0596   collectedStrips_.clear();
0597   collectedPixels_.clear();
0598 }
0599 
0600 #include "FWCore/PluginManager/interface/ModuleDef.h"
0601 #include "FWCore/Framework/interface/MakerMacros.h"
0602 DEFINE_FWK_MODULE(HITrackClusterRemover);