Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-01-22 07:34:25

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 "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0008 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0009 
0010 #include "DataFormats/SiStripCluster/interface/SiStripCluster.h"
0011 #include "DataFormats/SiPixelCluster/interface/SiPixelCluster.h"
0012 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2D.h"
0013 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit1D.h"
0014 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2D.h"
0015 #include "DataFormats/TrackerRecHit2D/interface/ProjectedSiStripRecHit2D.h"
0016 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
0017 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHitCollection.h"
0018 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2DCollection.h"
0019 #include "DataFormats/Common/interface/Handle.h"
0020 #include "DataFormats/Common/interface/DetSetVector.h"
0021 #include "DataFormats/Common/interface/DetSetVectorNew.h"
0022 #include "DataFormats/Common/interface/ContainerMask.h"
0023 #include "DataFormats/Provenance/interface/ProductID.h"
0024 
0025 #include "DataFormats/TrackReco/interface/Track.h"
0026 #include "DataFormats/TrackerRecHit2D/interface/ClusterRemovalInfo.h"
0027 
0028 #include "TrackingTools/PatternTools/interface/Trajectory.h"
0029 #include "TrackingTools/PatternTools/interface/TrajTrackAssociation.h"
0030 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
0031 
0032 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0033 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0034 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
0035 #include "Geometry/CommonDetUnit/interface/GeomDetType.h"
0036 
0037 //
0038 // class decleration
0039 //
0040 
0041 class HLTTrackClusterRemoverNew final : public edm::stream::EDProducer<> {
0042 public:
0043   HLTTrackClusterRemoverNew(const edm::ParameterSet &iConfig);
0044   ~HLTTrackClusterRemoverNew() override;
0045   void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override;
0046   static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
0047 
0048 private:
0049   edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> const tTrackerGeom_;
0050   struct ParamBlock {
0051     ParamBlock() : isSet_(false), usesCharge_(false) {}
0052     ParamBlock(const edm::ParameterSet &iConfig)
0053         : isSet_(true),
0054           maxChi2_(iConfig.getParameter<double>("maxChi2")),
0055           maxCharge_(iConfig.getParameter<double>("maxCharge")),
0056           minGoodPixelCharge_(iConfig.getParameter<double>("minGoodPixelCharge")),
0057           minGoodStripCharge_(iConfig.getParameter<double>("minGoodStripCharge")),
0058           maxSize_(iConfig.getParameter<uint32_t>("maxSize")),
0059           usesCharge_(maxCharge_ > 0.),
0060           usesSize_(maxSize_ > 0.),
0061           cutOnPixelCharge_(minGoodPixelCharge_ > 0.),
0062           cutOnStripCharge_(minGoodStripCharge_ > 0.) {}
0063     bool isSet_;
0064     float maxChi2_, maxCharge_, minGoodPixelCharge_, minGoodStripCharge_;
0065     size_t maxSize_;
0066     bool usesCharge_, usesSize_, cutOnPixelCharge_, cutOnStripCharge_;
0067 
0068     static void fillPSetDescription(edm::ParameterSetDescription &desc) {
0069       desc.addOptional<double>("maxChi2");
0070       desc.add<double>("maxCharge", 0.);
0071       desc.add<double>("minGoodPixelCharge", 0.);
0072       desc.add<double>("minGoodStripCharge", 0.);
0073       desc.add<uint32_t>("maxSize", 0.);
0074     }
0075   };
0076   static const unsigned int NumberOfParamBlocks = 6;
0077 
0078   bool doTracks_;
0079   bool doStrip_, doPixel_;
0080   bool mergeOld_;
0081 
0082   typedef edm::ContainerMask<edmNew::DetSetVector<SiPixelCluster> > PixelMaskContainer;
0083   typedef edm::ContainerMask<edmNew::DetSetVector<SiStripCluster> > StripMaskContainer;
0084   edm::EDGetTokenT<edmNew::DetSetVector<SiPixelCluster> > pixelClusters_;
0085   edm::EDGetTokenT<edmNew::DetSetVector<SiStripCluster> > stripClusters_;
0086   edm::EDGetTokenT<PixelMaskContainer> oldPxlMaskToken_;
0087   edm::EDGetTokenT<StripMaskContainer> oldStrMaskToken_;
0088   edm::EDGetTokenT<std::vector<Trajectory> > trajectories_;
0089 
0090   ParamBlock pblocks_[NumberOfParamBlocks];
0091   void readPSet(const edm::ParameterSet &iConfig,
0092                 const std::string &name,
0093                 int id1 = -1,
0094                 int id2 = -1,
0095                 int id3 = -1,
0096                 int id4 = -1,
0097                 int id5 = -1,
0098                 int id6 = -1);
0099 
0100   std::vector<uint8_t> pixels, strips;                  // avoid unneed alloc/dealloc of this
0101   edm::ProductID pixelSourceProdID, stripSourceProdID;  // ProdIDs refs must point to (for consistency tests)
0102 
0103   inline void process(const TrackingRecHit *hit, float chi2, const TrackerGeometry *tg);
0104   inline void process(const OmniClusterRef &cluRef, uint32_t subdet);
0105 
0106   template <typename T>
0107   std::unique_ptr<edmNew::DetSetVector<T> > cleanup(const edmNew::DetSetVector<T> &oldClusters,
0108                                                     const std::vector<uint8_t> &isGood,
0109                                                     reco::ClusterRemovalInfo::Indices &refs,
0110                                                     const reco::ClusterRemovalInfo::Indices *oldRefs);
0111 
0112   // Carries in full removal info about a given det from oldRefs
0113   void mergeOld(reco::ClusterRemovalInfo::Indices &refs, const reco::ClusterRemovalInfo::Indices &oldRefs);
0114 
0115   bool makeProducts_;
0116   bool doStripChargeCheck_, doPixelChargeCheck_;
0117   std::vector<bool> collectedRegStrips_;
0118   std::vector<bool> collectedPixels_;
0119 };
0120 
0121 using namespace std;
0122 using namespace edm;
0123 using namespace reco;
0124 
0125 void HLTTrackClusterRemoverNew::readPSet(
0126     const edm::ParameterSet &iConfig, const std::string &name, int id1, int id2, int id3, int id4, int id5, int id6) {
0127   if (iConfig.exists(name)) {
0128     ParamBlock pblock(iConfig.getParameter<ParameterSet>(name));
0129     if (id1 == -1) {
0130       fill(pblocks_, pblocks_ + NumberOfParamBlocks, pblock);
0131     } else {
0132       pblocks_[id1] = pblock;
0133       if (id2 != -1)
0134         pblocks_[id2] = pblock;
0135       if (id3 != -1)
0136         pblocks_[id3] = pblock;
0137       if (id4 != -1)
0138         pblocks_[id4] = pblock;
0139       if (id5 != -1)
0140         pblocks_[id5] = pblock;
0141       if (id6 != -1)
0142         pblocks_[id6] = pblock;
0143     }
0144   }
0145 }
0146 
0147 HLTTrackClusterRemoverNew::HLTTrackClusterRemoverNew(const ParameterSet &iConfig)
0148     : tTrackerGeom_(esConsumes<TrackerGeometry, TrackerDigiGeometryRecord>()),
0149       doTracks_(!iConfig.getParameter<edm::InputTag>("trajectories").label().empty()),
0150       doStrip_(iConfig.getParameter<bool>("doStrip")),
0151       doPixel_(iConfig.getParameter<bool>("doPixel")),
0152       mergeOld_(false),
0153       makeProducts_(true),
0154       doStripChargeCheck_(iConfig.getParameter<bool>("doStripChargeCheck")),
0155       doPixelChargeCheck_(iConfig.getParameter<bool>("doPixelChargeCheck"))
0156 
0157 {
0158   oldPxlMaskToken_ = consumes<PixelMaskContainer>(iConfig.getParameter<InputTag>("oldClusterRemovalInfo"));
0159   oldStrMaskToken_ = consumes<StripMaskContainer>(iConfig.getParameter<InputTag>("oldClusterRemovalInfo"));
0160   if (not(iConfig.getParameter<InputTag>("oldClusterRemovalInfo") == edm::InputTag())) {
0161     mergeOld_ = true;
0162   }
0163 
0164   if ((doPixelChargeCheck_ && !doPixel_) || (doStripChargeCheck_ && !doStrip_))
0165     throw cms::Exception("Configuration Error")
0166         << "HLTTrackClusterRemoverNew: Charge check asked without cluster collection ";
0167   if (doPixelChargeCheck_)
0168     throw cms::Exception("Configuration Error")
0169         << "HLTTrackClusterRemoverNew: Pixel cluster charge check not yet implemented";
0170 
0171   fill(pblocks_, pblocks_ + NumberOfParamBlocks, ParamBlock());
0172   readPSet(iConfig, "Common", -1);
0173   if (doPixel_) {
0174     readPSet(iConfig, "Pixel", 0, 1);
0175     readPSet(iConfig, "PXB", 0);
0176     readPSet(iConfig, "PXE", 1);
0177   }
0178   if (doStrip_) {
0179     readPSet(iConfig, "Strip", 2, 3, 4, 5);
0180     readPSet(iConfig, "StripInner", 2, 3);
0181     readPSet(iConfig, "StripOuter", 4, 5);
0182     readPSet(iConfig, "TIB", 2);
0183     readPSet(iConfig, "TID", 3);
0184     readPSet(iConfig, "TOB", 4);
0185     readPSet(iConfig, "TEC", 5);
0186   }
0187 
0188   bool usingCharge = false;
0189   for (size_t i = 0; i < NumberOfParamBlocks; ++i) {
0190     if (!pblocks_[i].isSet_)
0191       throw cms::Exception("Configuration Error")
0192           << "HLTTrackClusterRemoverNew: Missing configuration for detector with subDetID = " << (i + 1);
0193     if (pblocks_[i].usesCharge_ && !usingCharge) {
0194       throw cms::Exception("Configuration Error")
0195           << "HLTTrackClusterRemoverNew: Configuration for subDetID = " << (i + 1)
0196           << " uses cluster charge, which is not enabled.";
0197     }
0198   }
0199 
0200   //    trajectories_ = consumes<vector<Trajectory> >(iConfig.getParameter<InputTag>("trajectories"));
0201   if (doTracks_)
0202     trajectories_ = consumes<vector<Trajectory> >(iConfig.getParameter<InputTag>("trajectories"));
0203   if (doPixel_)
0204     pixelClusters_ = consumes<edmNew::DetSetVector<SiPixelCluster> >(iConfig.getParameter<InputTag>("pixelClusters"));
0205   if (doStrip_)
0206     stripClusters_ = consumes<edmNew::DetSetVector<SiStripCluster> >(iConfig.getParameter<InputTag>("stripClusters"));
0207   if (mergeOld_) {
0208     oldPxlMaskToken_ = consumes<PixelMaskContainer>(iConfig.getParameter<InputTag>("oldClusterRemovalInfo"));
0209     oldStrMaskToken_ = consumes<StripMaskContainer>(iConfig.getParameter<InputTag>("oldClusterRemovalInfo"));
0210   }
0211 
0212   //produces<edmNew::DetSetVector<SiPixelClusterRefNew> >();
0213   //produces<edmNew::DetSetVector<SiStripRecHit1D::ClusterRegionalRef> >();
0214 
0215   produces<edm::ContainerMask<edmNew::DetSetVector<SiPixelCluster> > >();
0216   produces<edm::ContainerMask<edmNew::DetSetVector<SiStripCluster> > >();
0217 }
0218 
0219 HLTTrackClusterRemoverNew::~HLTTrackClusterRemoverNew() {}
0220 
0221 void HLTTrackClusterRemoverNew::mergeOld(ClusterRemovalInfo::Indices &refs,
0222                                          const ClusterRemovalInfo::Indices &oldRefs) {
0223   for (size_t i = 0, n = refs.size(); i < n; ++i) {
0224     refs[i] = oldRefs[refs[i]];
0225   }
0226 }
0227 
0228 template <typename T>
0229 std::unique_ptr<edmNew::DetSetVector<T> > HLTTrackClusterRemoverNew::cleanup(
0230     const edmNew::DetSetVector<T> &oldClusters,
0231     const std::vector<uint8_t> &isGood,
0232     reco::ClusterRemovalInfo::Indices &refs,
0233     const reco::ClusterRemovalInfo::Indices *oldRefs) {
0234   typedef typename edmNew::DetSetVector<T> DSV;
0235   typedef typename edmNew::DetSetVector<T>::FastFiller DSF;
0236   typedef typename edmNew::DetSet<T> DS;
0237   auto output = std::make_unique<DSV>();
0238   output->reserve(oldClusters.size(), oldClusters.dataSize());
0239 
0240   // cluster removal loop
0241   const T *firstOffset = &oldClusters.data().front();
0242   for (typename DSV::const_iterator itdet = oldClusters.begin(), enddet = oldClusters.end(); itdet != enddet; ++itdet) {
0243     DS oldDS = *itdet;
0244 
0245     if (oldDS.empty())
0246       continue;  // skip empty detsets
0247 
0248     uint32_t id = oldDS.detId();
0249     DSF outds(*output, id);
0250 
0251     for (typename DS::const_iterator it = oldDS.begin(), ed = oldDS.end(); it != ed; ++it) {
0252       uint32_t index = ((&*it) - firstOffset);
0253       if (isGood[index]) {
0254         outds.push_back(*it);
0255         refs.push_back(index);
0256         //std::cout << "HLTTrackClusterRemoverNew::cleanup " << typeid(T).name() << " reference " << index << " to " << (refs.size() - 1) << std::endl;
0257       }
0258     }
0259     if (outds.empty())
0260       outds.abort();  // not write in an empty DSV
0261   }
0262   //    std::cout<<"fraction: "<<fraction<<std::endl;
0263   if (oldRefs != nullptr)
0264     mergeOld(refs, *oldRefs);
0265   return output;
0266 }
0267 
0268 void HLTTrackClusterRemoverNew::process(OmniClusterRef const &clusterReg, uint32_t subdet) {
0269   if (clusterReg.id() != stripSourceProdID)
0270     throw cms::Exception("Inconsistent Data")
0271         << "HLTTrackClusterRemoverNew: strip cluster ref from Product ID = " << clusterReg.id()
0272         << " does not match with source cluster collection (ID = " << stripSourceProdID << ")\n.";
0273 
0274   if (collectedRegStrips_.size() <= clusterReg.key()) {
0275     edm::LogError("BadCollectionSize") << collectedRegStrips_.size() << " is smaller than " << clusterReg.key();
0276 
0277     assert(collectedRegStrips_.size() > clusterReg.key());
0278   }
0279   collectedRegStrips_[clusterReg.key()] = true;
0280 }
0281 
0282 void HLTTrackClusterRemoverNew::process(const TrackingRecHit *hit, float chi2, const TrackerGeometry *tg) {
0283   DetId detid = hit->geographicalId();
0284   uint32_t subdet = detid.subdetId();
0285 
0286   assert((subdet > 0) && (subdet <= NumberOfParamBlocks));
0287 
0288   // chi2 cut
0289   if (chi2 > pblocks_[subdet - 1].maxChi2_)
0290     return;
0291 
0292   if (GeomDetEnumerators::isTrackerPixel(tg->geomDetSubDetector(subdet))) {
0293     //      std::cout<<"process pxl hit"<<std::endl;
0294     if (!doPixel_)
0295       return;
0296     // this is a pixel, and i *know* it is
0297     const SiPixelRecHit *pixelHit = static_cast<const SiPixelRecHit *>(hit);
0298 
0299     SiPixelRecHit::ClusterRef cluster = pixelHit->cluster();
0300     if (cluster.id() != pixelSourceProdID)
0301       throw cms::Exception("Inconsistent Data")
0302           << "HLTTrackClusterRemoverNew: pixel cluster ref from Product ID = " << cluster.id()
0303           << " does not match with source cluster collection (ID = " << pixelSourceProdID << ")\n.";
0304 
0305     assert(cluster.id() == pixelSourceProdID);
0306     //DBG// cout << "HIT NEW PIXEL DETID = " << detid.rawId() << ", Cluster [ " << cluster.key().first << " / " <<  cluster.key().second << " ] " << endl;
0307 
0308     // if requested, cut on cluster size
0309     if (pblocks_[subdet - 1].usesSize_ && (cluster->pixels().size() > pblocks_[subdet - 1].maxSize_))
0310       return;
0311 
0312     // mark as used
0313     //pixels[cluster.key()] = false;
0314     assert(collectedPixels_.size() > cluster.key());
0315     collectedPixels_[cluster.key()] = true;
0316   } else {  // aka Strip
0317     if (!doStrip_)
0318       return;
0319     const type_info &hitType = typeid(*hit);
0320     if (hitType == typeid(SiStripRecHit2D)) {
0321       const SiStripRecHit2D *stripHit = static_cast<const SiStripRecHit2D *>(hit);
0322       //DBG//     cout << "Plain RecHit 2D: " << endl;
0323       process(stripHit->omniClusterRef(), subdet);
0324       //      int clusCharge=0;
0325       //      for ( auto cAmp : stripHit->omniClusterRef().stripCluster().amplitudes() ) clusCharge+=cAmp;
0326       //      std::cout << "[HLTTrackClusterRemoverNew::process (SiStripRecHit2D) chi2: " << chi2 << " [" << subdet << " --> charge: " << clusCharge << "]" << std::endl;
0327     } else if (hitType == typeid(SiStripRecHit1D)) {
0328       const SiStripRecHit1D *hit1D = static_cast<const SiStripRecHit1D *>(hit);
0329       process(hit1D->omniClusterRef(), subdet);
0330       //      int clusCharge=0;
0331       //      for ( auto cAmp : hit1D->omniClusterRef().stripCluster().amplitudes() ) clusCharge+=cAmp;
0332       //      std::cout << "[HLTTrackClusterRemoverNew::process (SiStripRecHit1D) chi2: " << chi2 << " [" << subdet << " --> charge: " << clusCharge << "]" << std::endl;
0333     } else if (hitType == typeid(SiStripMatchedRecHit2D)) {
0334       const SiStripMatchedRecHit2D *matchHit = static_cast<const SiStripMatchedRecHit2D *>(hit);
0335       //DBG//     cout << "Matched RecHit 2D: " << endl;
0336       process(matchHit->monoClusterRef(), subdet);
0337       //        int clusCharge=0;
0338       //        for ( auto cAmp : matchHit->monoClusterRef().stripCluster().amplitudes() ) clusCharge+=cAmp;
0339       //        std::cout << "[HLTTrackClusterRemoverNew::process (SiStripMatchedRecHit2D:mono) chi2: " << chi2 << " [" << subdet << " --> charge: " << clusCharge << "]" << std::endl;
0340 
0341       process(matchHit->stereoClusterRef(), subdet);
0342       //        clusCharge=0;
0343       //        for ( auto cAmp : matchHit->stereoClusterRef().stripCluster().amplitudes() ) clusCharge+=cAmp;
0344       //        std::cout << "[HLTTrackClusterRemoverNew::process (SiStripMatchedRecHit2D:stereo) chi2: " << chi2 << " [" << subdet << " --> charge: " << clusCharge << "]" << std::endl;
0345 
0346     } else if (hitType == typeid(ProjectedSiStripRecHit2D)) {
0347       const ProjectedSiStripRecHit2D *projHit = static_cast<const ProjectedSiStripRecHit2D *>(hit);
0348       //DBG//     cout << "Projected RecHit 2D: " << endl;
0349       process(projHit->originalHit().omniClusterRef(), subdet);
0350       //        int clusCharge=0;
0351       //        for ( auto cAmp : projHit->originalHit().omniClusterRef().stripCluster().amplitudes() ) clusCharge+=cAmp;
0352       //        std::cout << "[HLTTrackClusterRemoverNew::process (ProjectedSiStripRecHit2D) chi2: " << chi2 << " [" << subdet << " --> charge: " << clusCharge << "]" << std::endl;
0353     } else
0354       throw cms::Exception("NOT IMPLEMENTED")
0355           << "Don't know how to handle " << hitType.name() << " on detid " << detid.rawId() << "\n";
0356   }
0357 }
0358 
0359 /*   Schematic picture of n-th step Iterative removal
0360  *   (that os removing clusters after n-th step tracking)
0361  *   clusters:    [ C1 ] -> [ C2 ] -> ... -> [ Cn ] -> [ Cn + 1 ]
0362  *                   ^                          ^           ^--- OUTPUT "new" ID
0363  *                   |-- before any removal     |----- Source clusters                                                                             
0364  *                   |-- OUTPUT "old" ID        |----- Hits in Traj. point here                       
0365  *                   |                          \----- Old ClusterRemovalInfo "new" ID                 
0366  *                   \-- Old ClusterRemovalInfo "old" ID                                                  
0367  */
0368 
0369 void HLTTrackClusterRemoverNew::produce(Event &iEvent, const EventSetup &iSetup) {
0370   ProductID pixelOldProdID, stripOldProdID;
0371 
0372   const auto &tgh = &iSetup.getData(tTrackerGeom_);
0373 
0374   edm::Handle<edmNew::DetSetVector<SiPixelCluster> > pixelClusters;
0375   if (doPixel_) {
0376     iEvent.getByToken(pixelClusters_, pixelClusters);
0377     pixelSourceProdID = pixelClusters.id();
0378   }
0379 
0380   edm::Handle<edmNew::DetSetVector<SiStripCluster> > stripClusters;
0381   if (doStrip_) {
0382     iEvent.getByToken(stripClusters_, stripClusters);
0383     stripSourceProdID = stripClusters.id();
0384   }
0385 
0386   //Handle<TrajTrackAssociationCollection> trajectories;
0387   edm::Handle<vector<Trajectory> > trajectories;
0388   iEvent.getByToken(trajectories_, trajectories);
0389 
0390   if (mergeOld_) {
0391     edm::Handle<PixelMaskContainer> oldPxlMask;
0392     edm::Handle<StripMaskContainer> oldStrMask;
0393     iEvent.getByToken(oldPxlMaskToken_, oldPxlMask);
0394     iEvent.getByToken(oldStrMaskToken_, oldStrMask);
0395     LogDebug("TrackClusterRemover") << "to merge in, " << oldStrMask->size() << " strp and " << oldPxlMask->size()
0396                                     << " pxl";
0397     oldStrMask->copyMaskTo(collectedRegStrips_);
0398     oldPxlMask->copyMaskTo(collectedPixels_);
0399     collectedRegStrips_.resize(stripClusters->dataSize(), false);
0400   } else {
0401     collectedRegStrips_.resize(stripClusters->dataSize(), false);
0402     collectedPixels_.resize(pixelClusters->dataSize(), false);
0403   }
0404 
0405   //for (TrajTrackAssociationCollection::const_iterator it = trajectories->begin(), ed = trajectories->end(); it != ed; ++it)  {
0406   //    const Trajectory &tj = * it->key;
0407 
0408   for (std::vector<Trajectory>::const_iterator it = trajectories->begin(), ed = trajectories->end(); it != ed; ++it) {
0409     const Trajectory &tj = *it;
0410     const std::vector<TrajectoryMeasurement> &tms = tj.measurements();
0411 
0412     std::vector<TrajectoryMeasurement>::const_iterator itm, endtm;
0413     for (itm = tms.begin(), endtm = tms.end(); itm != endtm; ++itm) {
0414       const TrackingRecHit *hit = itm->recHit()->hit();
0415       if (!hit->isValid())
0416         continue;
0417       //        std::cout<<"process hit"<<std::endl;
0418       process(hit, itm->estimate(), tgh);
0419     }
0420   }
0421 
0422   //    std::cout << " => collectedRegStrips_: " << collectedRegStrips_.size() << std::endl;
0423   //    std::cout << " total strip to skip (before charge check): "<<std::count(collectedRegStrips_.begin(),collectedRegStrips_.end(),true) << std::endl;
0424 
0425   // checks only cluster already found! creates regression
0426   if (doStripChargeCheck_) {
0427     //      std::cout << "[HLTTrackClusterRemoverNew::produce] doStripChargeCheck_: " << (doStripChargeCheck_ ? "true" : "false") << " stripClusters: " << stripClusters->size() << std::endl;
0428 
0429     auto const &clusters = stripClusters->data();
0430     for (auto const &item : stripClusters->ids()) {
0431       if (!item.isValid())
0432         continue;  // not umpacked
0433 
0434       DetId detid = item.id;
0435       uint32_t subdet = detid.subdetId();
0436       if (!pblocks_[subdet - 1].cutOnStripCharge_)
0437         continue;
0438 
0439       //    std::cout << " i: " << i << " --> detid: " << detid << " --> subdet: " << subdet << std::endl;
0440 
0441       for (int i = item.offset; i < item.offset + int(item.size); ++i) {
0442         int clusCharge = 0;
0443         for (auto cAmp : clusters[i].amplitudes())
0444           clusCharge += cAmp;
0445 
0446         //  if (clusCharge < pblocks_[subdet-1].minGoodStripCharge_) std::cout << " clusCharge: " << clusCharge << std::endl;
0447         if (clusCharge < pblocks_[subdet - 1].minGoodStripCharge_)
0448           collectedRegStrips_[i] = true;  // (|= does not work!)
0449       }
0450     }
0451   }
0452 
0453   //    std::cout << " => collectedRegStrips_: " << collectedRegStrips_.size() << std::endl;
0454   //    std::cout << " total strip to skip: "<<std::count(collectedRegStrips_.begin(),collectedRegStrips_.end(),true) << std::endl;
0455   LogDebug("TrackClusterRemover") << "total strip to skip: "
0456                                   << std::count(collectedRegStrips_.begin(), collectedRegStrips_.end(), true);
0457   // std::cout << "TrackClusterRemover" <<" total strip to skip: "<<std::count(collectedRegStrips_.begin(),collectedRegStrips_.end(),true)<<std::endl;
0458   iEvent.put(std::make_unique<StripMaskContainer>(edm::RefProd<edmNew::DetSetVector<SiStripCluster> >(stripClusters),
0459                                                   collectedRegStrips_));
0460 
0461   LogDebug("TrackClusterRemover") << "total pxl to skip: "
0462                                   << std::count(collectedPixels_.begin(), collectedPixels_.end(), true);
0463   iEvent.put(std::make_unique<PixelMaskContainer>(edm::RefProd<edmNew::DetSetVector<SiPixelCluster> >(pixelClusters),
0464                                                   collectedPixels_));
0465 
0466   collectedRegStrips_.clear();
0467   collectedPixels_.clear();
0468 }
0469 
0470 void HLTTrackClusterRemoverNew::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
0471   edm::ParameterSetDescription desc;
0472 
0473   // Define the structure for optional ParameterSets like Common, Pixel, Strip, etc.
0474   edm::ParameterSetDescription paramBlockDesc;
0475   ParamBlock::fillPSetDescription(paramBlockDesc);
0476 
0477   // Add all possible ParameterSets (optional) using the above structure
0478   desc.addOptional<edm::ParameterSetDescription>("Common", paramBlockDesc);
0479   desc.addOptional<edm::ParameterSetDescription>("Pixel", paramBlockDesc);
0480   desc.addOptional<edm::ParameterSetDescription>("Strip", paramBlockDesc);
0481   desc.addOptional<edm::ParameterSetDescription>("PXB", paramBlockDesc);
0482   desc.addOptional<edm::ParameterSetDescription>("PXE", paramBlockDesc);
0483   desc.addOptional<edm::ParameterSetDescription>("StripInner", paramBlockDesc);
0484   desc.addOptional<edm::ParameterSetDescription>("StripOuter", paramBlockDesc);
0485   desc.addOptional<edm::ParameterSetDescription>("TIB", paramBlockDesc);
0486   desc.addOptional<edm::ParameterSetDescription>("TID", paramBlockDesc);
0487   desc.addOptional<edm::ParameterSetDescription>("TOB", paramBlockDesc);
0488   desc.addOptional<edm::ParameterSetDescription>("TEC", paramBlockDesc);
0489 
0490   desc.add<bool>("doStrip", true);
0491   desc.add<bool>("doPixel", true);
0492   desc.add<bool>("doStripChargeCheck", false);
0493   desc.add<bool>("doPixelChargeCheck", false);
0494   desc.add<InputTag>("trajectories", edm::InputTag(""));
0495   desc.add<InputTag>("pixelClusters");
0496   desc.add<InputTag>("stripClusters");
0497   desc.add<InputTag>("oldClusterRemovalInfo", edm::InputTag(""));
0498   descriptions.addDefault(desc);
0499 }
0500 
0501 #include "FWCore/PluginManager/interface/ModuleDef.h"
0502 #include "FWCore/Framework/interface/MakerMacros.h"
0503 DEFINE_FWK_MODULE(HLTTrackClusterRemoverNew);