File indexing completed on 2024-04-06 12:26:33
0001 #include "FWCore/Framework/interface/global/EDProducer.h"
0002 #include "FWCore/Framework/interface/Event.h"
0003 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0004 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0005 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0006
0007 #include "FWCore/Utilities/interface/InputTag.h"
0008
0009 #include "DataFormats/SiPixelCluster/interface/SiPixelCluster.h"
0010 #include "DataFormats/Phase2TrackerCluster/interface/Phase2TrackerCluster1D.h"
0011
0012 #include "DataFormats/Common/interface/ValueMap.h"
0013 #include "DataFormats/Common/interface/DetSetVectorNew.h"
0014 #include "DataFormats/Provenance/interface/ProductID.h"
0015 #include "DataFormats/Common/interface/ContainerMask.h"
0016
0017 #include "DataFormats/DetId/interface/DetId.h"
0018
0019 #include "DataFormats/TrackReco/interface/Track.h"
0020 #include "DataFormats/TrackerRecHit2D/interface/ClusterRemovalInfo.h"
0021 #include "DataFormats/TrackerRecHit2D/interface/VectorHit.h"
0022
0023 #include "TrackingTools/PatternTools/interface/TrackCollectionTokens.h"
0024
0025 #include "RecoTracker/TransientTrackingRecHit/interface/Traj2TrackHits.h"
0026
0027 #include <limits>
0028
0029
0030
0031
0032
0033
0034
0035 namespace {
0036
0037 class TrackClusterRemoverPhase2 final : public edm::global::EDProducer<> {
0038 public:
0039 TrackClusterRemoverPhase2(const edm::ParameterSet& iConfig);
0040 ~TrackClusterRemoverPhase2() override {}
0041 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0042
0043 private:
0044 void produce(edm::StreamID, edm::Event& evt, const edm::EventSetup&) const override;
0045
0046 using PixelMaskContainer = edm::ContainerMask<edmNew::DetSetVector<SiPixelCluster>>;
0047 using Phase2OTMaskContainer = edm::ContainerMask<edmNew::DetSetVector<Phase2TrackerCluster1D>>;
0048
0049 using QualityMaskCollection = std::vector<unsigned char>;
0050
0051 const unsigned char maxChi2_;
0052 const int minNumberOfLayersWithMeasBeforeFiltering_;
0053 const reco::TrackBase::TrackQuality trackQuality_;
0054
0055 const TrackCollectionTokens tracks_;
0056 edm::EDGetTokenT<QualityMaskCollection> srcQuals;
0057
0058 const edm::EDGetTokenT<edmNew::DetSetVector<SiPixelCluster>> pixelClusters_;
0059 const edm::EDGetTokenT<edmNew::DetSetVector<Phase2TrackerCluster1D>> phase2OTClusters_;
0060
0061 edm::EDGetTokenT<PixelMaskContainer> oldPxlMaskToken_;
0062 edm::EDGetTokenT<Phase2OTMaskContainer> oldPh2OTMaskToken_;
0063
0064
0065 edm::EDGetTokenT<edm::ValueMap<int>> overrideTrkQuals_;
0066 };
0067
0068 void TrackClusterRemoverPhase2::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0069 edm::ParameterSetDescription desc;
0070 desc.add<edm::InputTag>("trajectories", edm::InputTag());
0071 desc.add<edm::InputTag>("trackClassifier", edm::InputTag("", "QualityMasks"));
0072 desc.add<edm::InputTag>("phase2pixelClusters", edm::InputTag("siPixelClusters"));
0073 desc.add<edm::InputTag>("phase2OTClusters", edm::InputTag("siPhase2Clusters"));
0074 desc.add<edm::InputTag>("oldClusterRemovalInfo", edm::InputTag());
0075
0076 desc.add<std::string>("TrackQuality", "highPurity");
0077 desc.add<double>("maxChi2", 30.);
0078 desc.add<int>("minNumberOfLayersWithMeasBeforeFiltering", 0);
0079
0080 desc.add<edm::InputTag>("overrideTrkQuals", edm::InputTag());
0081
0082 descriptions.add("phase2trackClusterRemover", desc);
0083 }
0084
0085 TrackClusterRemoverPhase2::TrackClusterRemoverPhase2(const edm::ParameterSet& iConfig)
0086 : maxChi2_(Traj2TrackHits::toChi2x5(iConfig.getParameter<double>("maxChi2"))),
0087 minNumberOfLayersWithMeasBeforeFiltering_(
0088 iConfig.getParameter<int>("minNumberOfLayersWithMeasBeforeFiltering")),
0089 trackQuality_(reco::TrackBase::qualityByName(iConfig.getParameter<std::string>("TrackQuality"))),
0090
0091 tracks_(iConfig.getParameter<edm::InputTag>("trajectories"), consumesCollector()),
0092 pixelClusters_(
0093 consumes<edmNew::DetSetVector<SiPixelCluster>>(iConfig.getParameter<edm::InputTag>("phase2pixelClusters"))),
0094 phase2OTClusters_(consumes<edmNew::DetSetVector<Phase2TrackerCluster1D>>(
0095 iConfig.getParameter<edm::InputTag>("phase2OTClusters"))) {
0096 produces<edm::ContainerMask<edmNew::DetSetVector<SiPixelCluster>>>();
0097 produces<edm::ContainerMask<edmNew::DetSetVector<Phase2TrackerCluster1D>>>();
0098
0099
0100 auto const& overrideTrkQuals = iConfig.getParameter<edm::InputTag>("overrideTrkQuals");
0101 if (!overrideTrkQuals.label().empty())
0102 overrideTrkQuals_ = consumes<edm::ValueMap<int>>(overrideTrkQuals);
0103
0104 auto const& classifier = iConfig.getParameter<edm::InputTag>("trackClassifier");
0105 if (!classifier.label().empty())
0106 srcQuals = consumes<QualityMaskCollection>(classifier);
0107
0108 auto const& oldClusterRemovalInfo = iConfig.getParameter<edm::InputTag>("oldClusterRemovalInfo");
0109 if (!oldClusterRemovalInfo.label().empty()) {
0110 oldPxlMaskToken_ = consumes<PixelMaskContainer>(oldClusterRemovalInfo);
0111 oldPh2OTMaskToken_ = consumes<Phase2OTMaskContainer>(oldClusterRemovalInfo);
0112 }
0113 }
0114
0115 void TrackClusterRemoverPhase2::produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup&) const {
0116 edm::Handle<edmNew::DetSetVector<SiPixelCluster>> pixelClusters;
0117 iEvent.getByToken(pixelClusters_, pixelClusters);
0118 edm::Handle<edmNew::DetSetVector<Phase2TrackerCluster1D>> phase2OTClusters;
0119 iEvent.getByToken(phase2OTClusters_, phase2OTClusters);
0120
0121 std::vector<bool> collectedPixels;
0122 std::vector<bool> collectedPhase2OTs;
0123
0124 if (!oldPxlMaskToken_.isUninitialized()) {
0125 edm::Handle<PixelMaskContainer> oldPxlMask;
0126 edm::Handle<Phase2OTMaskContainer> oldPh2OTMask;
0127 iEvent.getByToken(oldPxlMaskToken_, oldPxlMask);
0128 iEvent.getByToken(oldPh2OTMaskToken_, oldPh2OTMask);
0129 LogDebug("TrackClusterRemoverPhase2")
0130 << "to merge in, " << oldPxlMask->size() << " phase2 pixel and " << oldPh2OTMask->size() << " phase2 OT";
0131 oldPxlMask->copyMaskTo(collectedPixels);
0132 oldPh2OTMask->copyMaskTo(collectedPhase2OTs);
0133 assert(phase2OTClusters->dataSize() >= collectedPhase2OTs.size());
0134 collectedPhase2OTs.resize(phase2OTClusters->dataSize(), false);
0135
0136 } else {
0137 collectedPixels.resize(pixelClusters->dataSize(), false);
0138 collectedPhase2OTs.resize(phase2OTClusters->dataSize(), false);
0139 }
0140
0141
0142
0143 unsigned char qualMask = ~0;
0144 if (trackQuality_ != reco::TrackBase::undefQuality)
0145 qualMask = 1 << trackQuality_;
0146
0147 auto const& tracks = tracks_.tracks(iEvent);
0148 auto s = tracks.size();
0149
0150 QualityMaskCollection oldStyle;
0151 QualityMaskCollection const* pquals = nullptr;
0152
0153 if (!overrideTrkQuals_.isUninitialized()) {
0154 edm::Handle<edm::ValueMap<int>> quals;
0155 iEvent.getByToken(overrideTrkQuals_, quals);
0156 assert(s == (*quals).size());
0157
0158 oldStyle.resize(s, 0);
0159 for (auto i = 0U; i < s; ++i)
0160 if ((*quals).get(i) > 0)
0161 oldStyle[i] = (255) & (*quals).get(i);
0162 pquals = &oldStyle;
0163 }
0164
0165 if (!srcQuals.isUninitialized()) {
0166 edm::Handle<QualityMaskCollection> hqual;
0167 iEvent.getByToken(srcQuals, hqual);
0168 pquals = hqual.product();
0169 }
0170
0171 for (auto i = 0U; i < s; ++i) {
0172 const reco::Track& track = tracks[i];
0173 bool goodTk = (pquals) ? (*pquals)[i] & qualMask : track.quality(trackQuality_);
0174 if (!goodTk)
0175 continue;
0176 if (track.hitPattern().trackerLayersWithMeasurement() < minNumberOfLayersWithMeasBeforeFiltering_)
0177 continue;
0178 auto const& chi2sX5 = track.extra()->chi2sX5();
0179 assert(chi2sX5.size() == track.recHitsSize());
0180 auto hb = track.recHitsBegin();
0181 for (unsigned int h = 0; h < track.recHitsSize(); h++) {
0182 auto const hit = *(hb + h);
0183 if (!hit->isValid())
0184 continue;
0185 if (chi2sX5[h] > maxChi2_)
0186 continue;
0187 auto const& thit = static_cast<BaseTrackerRecHit const&>(*hit);
0188 auto const& cluster = thit.firstClusterRef();
0189
0190 if (cluster.isPixel())
0191 collectedPixels[cluster.key()] = true;
0192 else if (cluster.isPhase2())
0193 collectedPhase2OTs[cluster.key()] = true;
0194
0195
0196 auto pVectorHits = dynamic_cast<VectorHit const*>(hit);
0197 if (pVectorHits != nullptr) {
0198 auto const& vectorHit = *pVectorHits;
0199 auto const& lowCluster = vectorHit.lowerClusterRef();
0200 auto const& uppCluster = vectorHit.upperClusterRef();
0201 LogTrace("TrackClusterRemoverPhase2")
0202 << "masking a VHit with lowCluster key: " << lowCluster.key() << " and upper key: " << uppCluster.key();
0203 if (lowCluster.isPhase2())
0204 collectedPhase2OTs[lowCluster.key()] = true;
0205 if (uppCluster.isPhase2())
0206 collectedPhase2OTs[uppCluster.key()] = true;
0207 } else {
0208 LogTrace("TrackClusterRemoverPhase2") << "it is not a VHit.";
0209 }
0210 }
0211 }
0212
0213 auto removedPixelClusterMask = std::make_unique<PixelMaskContainer>(
0214 edm::RefProd<edmNew::DetSetVector<SiPixelCluster>>(pixelClusters), collectedPixels);
0215 LogDebug("TrackClusterRemoverPhase2")
0216 << "total pxl to skip: " << std::count(collectedPixels.begin(), collectedPixels.end(), true);
0217 iEvent.put(std::move(removedPixelClusterMask));
0218
0219 auto removedPhase2OTClusterMask = std::make_unique<Phase2OTMaskContainer>(
0220 edm::RefProd<edmNew::DetSetVector<Phase2TrackerCluster1D>>(phase2OTClusters), collectedPhase2OTs);
0221 LogDebug("TrackClusterRemoverPhase2")
0222 << "total ph2OT to skip: " << std::count(collectedPhase2OTs.begin(), collectedPhase2OTs.end(), true);
0223 iEvent.put(std::move(removedPhase2OTClusterMask));
0224 }
0225
0226 }
0227
0228 #include "FWCore/PluginManager/interface/ModuleDef.h"
0229 #include "FWCore/Framework/interface/MakerMacros.h"
0230 DEFINE_FWK_MODULE(TrackClusterRemoverPhase2);