File indexing completed on 2024-04-06 12:11:23
0001
0002 #include <memory>
0003
0004
0005 #include "FWCore/Framework/interface/Frameworkfwd.h"
0006 #include "FWCore/Framework/interface/stream/EDProducer.h"
0007 #include "FWCore/Framework/interface/Event.h"
0008 #include "FWCore/Framework/interface/MakerMacros.h"
0009 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0010 #include "FWCore/Framework/interface/ESHandle.h"
0011
0012
0013 #include "DataFormats/Common/interface/ValueMap.h"
0014 #include "DataFormats/TrackReco/interface/Track.h"
0015 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0016 #include "DataFormats/TrackerRecHit2D/interface/FastTrackerRecHit.h"
0017 #include "DataFormats/TrackerRecHit2D/interface/FastTrackerRecHitCollection.h"
0018
0019
0020 #include "TrackingTools/PatternTools/interface/TrackCollectionTokens.h"
0021
0022 class FastTrackerRecHitMaskProducer : public edm::stream::EDProducer<> {
0023 public:
0024 explicit FastTrackerRecHitMaskProducer(const edm::ParameterSet&);
0025
0026 ~FastTrackerRecHitMaskProducer() override {}
0027
0028 void produce(edm::Event&, const edm::EventSetup&) override;
0029
0030 private:
0031
0032 using QualityMaskCollection = std::vector<unsigned char>;
0033
0034
0035 const int minNumberOfLayersWithMeasBeforeFiltering_;
0036 const reco::TrackBase::TrackQuality trackQuality_;
0037
0038 const TrackCollectionTokens trajectories_;
0039 edm::EDGetTokenT<QualityMaskCollection> srcQuals;
0040
0041 edm::EDGetTokenT<FastTrackerRecHitCollection> recHits_;
0042
0043 edm::EDGetTokenT<std::vector<bool> > oldHitMaskToken_;
0044 };
0045
0046 FastTrackerRecHitMaskProducer::FastTrackerRecHitMaskProducer(const edm::ParameterSet& iConfig)
0047 : minNumberOfLayersWithMeasBeforeFiltering_(iConfig.getParameter<int>("minNumberOfLayersWithMeasBeforeFiltering")),
0048 trackQuality_(reco::TrackBase::qualityByName(iConfig.getParameter<std::string>("TrackQuality"))),
0049 trajectories_(iConfig.getParameter<edm::InputTag>("trajectories"), consumesCollector()),
0050 recHits_(consumes<FastTrackerRecHitCollection>(iConfig.getParameter<edm::InputTag>("recHits"))) {
0051 produces<std::vector<bool> >();
0052
0053 auto const& classifier = iConfig.getParameter<edm::InputTag>("trackClassifier");
0054 if (!classifier.label().empty())
0055 srcQuals = consumes<QualityMaskCollection>(classifier);
0056
0057 auto const& oldHitRemovalInfo = iConfig.getParameter<edm::InputTag>("oldHitRemovalInfo");
0058 if (!oldHitRemovalInfo.label().empty()) {
0059 oldHitMaskToken_ = consumes<std::vector<bool> >(oldHitRemovalInfo);
0060 }
0061 }
0062
0063 void FastTrackerRecHitMaskProducer::produce(edm::Event& iEvent, const edm::EventSetup& es) {
0064
0065
0066 std::unique_ptr<std::vector<bool> > collectedHits(new std::vector<bool>());
0067
0068
0069
0070 edm::Handle<FastTrackerRecHitCollection> recHits;
0071 iEvent.getByToken(recHits_, recHits);
0072
0073 if (!oldHitMaskToken_.isUninitialized()) {
0074 edm::Handle<std::vector<bool> > oldHitMasks;
0075 iEvent.getByToken(oldHitMaskToken_, oldHitMasks);
0076 collectedHits->insert(collectedHits->begin(), oldHitMasks->begin(), oldHitMasks->end());
0077 }
0078 collectedHits->resize(recHits->size(), false);
0079
0080 auto const& tracks = trajectories_.tracks(iEvent);
0081
0082 QualityMaskCollection const* pquals = nullptr;
0083 if (!srcQuals.isUninitialized()) {
0084 edm::Handle<QualityMaskCollection> hqual;
0085 iEvent.getByToken(srcQuals, hqual);
0086 pquals = hqual.product();
0087 }
0088
0089
0090 unsigned char qualMask = ~0;
0091 if (trackQuality_ != reco::TrackBase::undefQuality)
0092 qualMask = 1 << trackQuality_;
0093
0094
0095 for (auto i = 0U; i < tracks.size(); ++i) {
0096 const reco::Track& track = tracks[i];
0097 bool goodTk = (pquals) ? (*pquals)[i] & qualMask : track.quality(trackQuality_);
0098 if (!goodTk)
0099 continue;
0100 if (track.hitPattern().trackerLayersWithMeasurement() < minNumberOfLayersWithMeasBeforeFiltering_)
0101 continue;
0102 for (auto hitIt = track.recHitsBegin(); hitIt != track.recHitsEnd(); ++hitIt) {
0103 if (!(*hitIt)->isValid())
0104 continue;
0105 const FastTrackerRecHit& hit = static_cast<const FastTrackerRecHit&>(*(*hitIt));
0106
0107 for (unsigned id_index = 0; id_index < hit.nIds(); id_index++) {
0108 (*collectedHits)[unsigned(hit.id(id_index))] = true;
0109 }
0110 }
0111 }
0112
0113 iEvent.put(std::move(collectedHits));
0114 }
0115
0116 DEFINE_FWK_MODULE(FastTrackerRecHitMaskProducer);