CMS 3D CMS Logo

TrackClusterRemoverPhase2.cc
Go to the documentation of this file.
6 
8 
11 
16 
18 
22 
24 
26 
27 #include <limits>
28 
29 /* This is a copy of the TrackClusterRemover
30  * for Phase2 Tk upgrade
31  * FIXME:: changing with new phase2 pixel DataFormats!
32  * FIXME:: still to be factorized
33  */
34 
35 namespace {
36 
37  class TrackClusterRemoverPhase2 final : public edm::global::EDProducer<> {
38  public:
39  TrackClusterRemoverPhase2(const edm::ParameterSet& iConfig);
40  ~TrackClusterRemoverPhase2() override {}
41  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
42 
43  private:
44  void produce(edm::StreamID, edm::Event& evt, const edm::EventSetup&) const override;
45 
48 
49  using QualityMaskCollection = std::vector<unsigned char>;
50 
51  const unsigned char maxChi2_;
53  const reco::TrackBase::TrackQuality trackQuality_;
54 
55  const TrackCollectionTokens tracks_;
57 
60 
61  edm::EDGetTokenT<PixelMaskContainer> oldPxlMaskToken_;
62  edm::EDGetTokenT<Phase2OTMaskContainer> oldPh2OTMaskToken_;
63 
64  // backward compatibility during transition period
65  edm::EDGetTokenT<edm::ValueMap<int>> overrideTrkQuals_;
66  };
67 
70  desc.add<edm::InputTag>("trajectories", edm::InputTag());
71  desc.add<edm::InputTag>("trackClassifier", edm::InputTag("", "QualityMasks"));
72  desc.add<edm::InputTag>("phase2pixelClusters", edm::InputTag("siPixelClusters"));
73  desc.add<edm::InputTag>("phase2OTClusters", edm::InputTag("siPhase2Clusters"));
74  desc.add<edm::InputTag>("oldClusterRemovalInfo", edm::InputTag());
75 
76  desc.add<std::string>("TrackQuality", "highPurity");
77  desc.add<double>("maxChi2", 30.);
78  desc.add<int>("minNumberOfLayersWithMeasBeforeFiltering", 0);
79  // old mode
80  desc.add<edm::InputTag>("overrideTrkQuals", edm::InputTag());
81 
82  descriptions.add("phase2trackClusterRemover", desc);
83  }
84 
85  TrackClusterRemoverPhase2::TrackClusterRemoverPhase2(const edm::ParameterSet& iConfig)
86  : maxChi2_(Traj2TrackHits::toChi2x5(iConfig.getParameter<double>("maxChi2"))),
88  iConfig.getParameter<int>("minNumberOfLayersWithMeasBeforeFiltering")),
89  trackQuality_(reco::TrackBase::qualityByName(iConfig.getParameter<std::string>("TrackQuality"))),
90 
91  tracks_(iConfig.getParameter<edm::InputTag>("trajectories"), consumesCollector()),
92  pixelClusters_(
93  consumes<edmNew::DetSetVector<SiPixelCluster>>(iConfig.getParameter<edm::InputTag>("phase2pixelClusters"))),
94  phase2OTClusters_(consumes<edmNew::DetSetVector<Phase2TrackerCluster1D>>(
95  iConfig.getParameter<edm::InputTag>("phase2OTClusters"))) {
96  produces<edm::ContainerMask<edmNew::DetSetVector<SiPixelCluster>>>();
97  produces<edm::ContainerMask<edmNew::DetSetVector<Phase2TrackerCluster1D>>>();
98 
99  // old mode
100  auto const& overrideTrkQuals = iConfig.getParameter<edm::InputTag>("overrideTrkQuals");
101  if (!overrideTrkQuals.label().empty())
102  overrideTrkQuals_ = consumes<edm::ValueMap<int>>(overrideTrkQuals);
103 
104  auto const& classifier = iConfig.getParameter<edm::InputTag>("trackClassifier");
105  if (!classifier.label().empty())
106  srcQuals = consumes<QualityMaskCollection>(classifier);
107 
108  auto const& oldClusterRemovalInfo = iConfig.getParameter<edm::InputTag>("oldClusterRemovalInfo");
109  if (!oldClusterRemovalInfo.label().empty()) {
110  oldPxlMaskToken_ = consumes<PixelMaskContainer>(oldClusterRemovalInfo);
111  oldPh2OTMaskToken_ = consumes<Phase2OTMaskContainer>(oldClusterRemovalInfo);
112  }
113  }
114 
115  void TrackClusterRemoverPhase2::produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup&) const {
117  iEvent.getByToken(pixelClusters_, pixelClusters);
119  iEvent.getByToken(phase2OTClusters_, phase2OTClusters);
120 
121  std::vector<bool> collectedPixels;
122  std::vector<bool> collectedPhase2OTs;
123 
124  if (!oldPxlMaskToken_.isUninitialized()) {
127  iEvent.getByToken(oldPxlMaskToken_, oldPxlMask);
128  iEvent.getByToken(oldPh2OTMaskToken_, oldPh2OTMask);
129  LogDebug("TrackClusterRemoverPhase2")
130  << "to merge in, " << oldPxlMask->size() << " phase2 pixel and " << oldPh2OTMask->size() << " phase2 OT";
131  oldPxlMask->copyMaskTo(collectedPixels);
132  oldPh2OTMask->copyMaskTo(collectedPhase2OTs);
133  assert(phase2OTClusters->dataSize() >= collectedPhase2OTs.size());
134  collectedPhase2OTs.resize(phase2OTClusters->dataSize(), false);
135 
136  } else {
137  collectedPixels.resize(pixelClusters->dataSize(), false);
138  collectedPhase2OTs.resize(phase2OTClusters->dataSize(), false);
139  }
140 
141  // loop over trajectories, filter, mask clusters../
142 
143  unsigned char qualMask = ~0;
144  if (trackQuality_ != reco::TrackBase::undefQuality)
145  qualMask = 1 << trackQuality_;
146 
147  auto const& tracks = tracks_.tracks(iEvent);
148  auto s = tracks.size();
149 
150  QualityMaskCollection oldStyle;
151  QualityMaskCollection const* pquals = nullptr;
152 
153  if (!overrideTrkQuals_.isUninitialized()) {
155  iEvent.getByToken(overrideTrkQuals_, quals);
156  assert(s == (*quals).size());
157 
158  oldStyle.resize(s, 0);
159  for (auto i = 0U; i < s; ++i)
160  if ((*quals).get(i) > 0)
161  oldStyle[i] = (255) & (*quals).get(i);
162  pquals = &oldStyle;
163  }
164 
165  if (!srcQuals.isUninitialized()) {
167  iEvent.getByToken(srcQuals, hqual);
168  pquals = hqual.product();
169  }
170 
171  for (auto i = 0U; i < s; ++i) {
172  const reco::Track& track = tracks[i];
173  bool goodTk = (pquals) ? (*pquals)[i] & qualMask : track.quality(trackQuality_);
174  if (!goodTk)
175  continue;
176  if (track.hitPattern().trackerLayersWithMeasurement() < minNumberOfLayersWithMeasBeforeFiltering_)
177  continue;
178  auto const& chi2sX5 = track.extra()->chi2sX5();
179  assert(chi2sX5.size() == track.recHitsSize());
180  auto hb = track.recHitsBegin();
181  for (unsigned int h = 0; h < track.recHitsSize(); h++) {
182  auto const hit = *(hb + h);
183  if (!hit->isValid())
184  continue;
185  if (chi2sX5[h] > maxChi2_)
186  continue; // skip outliers
187  auto const& thit = static_cast<BaseTrackerRecHit const&>(*hit);
188  auto const& cluster = thit.firstClusterRef();
189  // FIXME when we will get also Phase2 pixel
190  if (cluster.isPixel())
191  collectedPixels[cluster.key()] = true;
192  else if (cluster.isPhase2())
193  collectedPhase2OTs[cluster.key()] = true;
194 
195  // Phase 2 OT is defined as Pixel detector (for now)
196  auto pVectorHits = dynamic_cast<VectorHit const*>(hit);
197  if (pVectorHits != nullptr) {
198  auto const& vectorHit = *pVectorHits;
199  auto const& lowCluster = vectorHit.lowerClusterRef();
200  auto const& uppCluster = vectorHit.upperClusterRef();
201  LogTrace("TrackClusterRemoverPhase2")
202  << "masking a VHit with lowCluster key: " << lowCluster.key() << " and upper key: " << uppCluster.key();
203  if (lowCluster.isPhase2())
204  collectedPhase2OTs[lowCluster.key()] = true;
205  if (uppCluster.isPhase2())
206  collectedPhase2OTs[uppCluster.key()] = true;
207  } else {
208  LogTrace("TrackClusterRemoverPhase2") << "it is not a VHit.";
209  }
210  }
211  }
212 
213  auto removedPixelClusterMask = std::make_unique<PixelMaskContainer>(
215  LogDebug("TrackClusterRemoverPhase2")
216  << "total pxl to skip: " << std::count(collectedPixels.begin(), collectedPixels.end(), true);
217  iEvent.put(std::move(removedPixelClusterMask));
218 
219  auto removedPhase2OTClusterMask = std::make_unique<Phase2OTMaskContainer>(
221  LogDebug("TrackClusterRemoverPhase2")
222  << "total ph2OT to skip: " << std::count(collectedPhase2OTs.begin(), collectedPhase2OTs.end(), true);
223  iEvent.put(std::move(removedPhase2OTClusterMask));
224  }
225 
226 } // namespace
227 
230 DEFINE_FWK_MODULE(TrackClusterRemoverPhase2);
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
Quality qualityByName(std::string const &name)
TrackQuality
track quality
Definition: TrackBase.h:150
T const * product() const
Definition: Handle.h:70
assert(be >=bs)
#define LogTrace(id)
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static void fillDescriptions(ConfigurationDescriptions &descriptions)
std::vector< unsigned char > QualityMaskCollection
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
virtual void produce(StreamID, Event &, EventSetup const &) const =0
Pixel cluster – collection of neighboring pixels above threshold.
fixed size matrix
HLT enums.
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
def move(src, dest)
Definition: eostools.py:511
#define LogDebug(id)