CMS 3D CMS Logo

TrackCollectionMerger.cc
Go to the documentation of this file.
2 
5 
12 
20 
23 
24 #include <vector>
25 #include <memory>
26 #include <cassert>
27 namespace {
28  class TrackCollectionMerger final : public edm::global::EDProducer<> {
29  public:
30  explicit TrackCollectionMerger(const edm::ParameterSet& conf)
31  : collectionCloner(producesCollector(), conf, true),
32  priorityName_(conf.getParameter<std::string>("trackAlgoPriorityOrder")),
33  m_foundHitBonus(conf.getParameter<double>("foundHitBonus")),
34  m_lostHitPenalty(conf.getParameter<double>("lostHitPenalty")),
35  m_shareFrac(conf.getParameter<double>("shareFrac")),
36  m_minShareHits(conf.getParameter<unsigned int>("minShareHits")),
37  m_minQuality(reco::TrackBase::qualityByName(conf.getParameter<std::string>("minQuality"))),
38  m_allowFirstHitShare(conf.getParameter<bool>("allowFirstHitShare")),
39  m_enableMerging(conf.getParameter<bool>("enableMerging")) {
40  for (auto const& it : conf.getParameter<std::vector<edm::InputTag>>("trackProducers"))
41  srcColls.emplace_back(it, consumesCollector());
42  for (auto const& it : conf.getParameter<std::vector<std::string>>("inputClassifiers")) {
43  srcMVAs.push_back(consumes<MVACollection>(edm::InputTag(it, "MVAValues")));
44  srcQuals.push_back(consumes<QualityMaskCollection>(edm::InputTag(it, "QualityMasks")));
45  }
46 
47  assert(srcColls.size() == srcQuals.size());
48 
49  produces<MVACollection>("MVAValues");
50  produces<QualityMaskCollection>("QualityMasks");
51  }
52 
53  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
55  desc.add<std::vector<edm::InputTag>>("trackProducers", std::vector<edm::InputTag>());
56  desc.add<std::vector<std::string>>("inputClassifiers", std::vector<std::string>());
57  desc.add<std::string>("trackAlgoPriorityOrder", "trackAlgoPriorityOrder");
58  desc.add<double>("shareFrac", .19);
59  desc.add<double>("foundHitBonus", 10.);
60  desc.add<double>("lostHitPenalty", 5.);
61  desc.add<unsigned int>("minShareHits", 2);
62  desc.add<bool>("allowFirstHitShare", true);
63  desc.add<bool>("enableMerging", true);
64  desc.add<std::string>("minQuality", "loose");
66  descriptions.add("TrackCollectionMerger", desc);
67  }
68 
69  private:
70  TrackCollectionCloner collectionCloner;
71  std::vector<TrackCollectionCloner::Tokens> srcColls;
72 
73  using MVACollection = std::vector<float>;
74  using QualityMaskCollection = std::vector<unsigned char>;
75 
76  using IHit = std::pair<unsigned int, TrackingRecHit const*>;
77  using IHitV = std::vector<IHit>;
78 
79  std::vector<edm::EDGetTokenT<MVACollection>> srcMVAs;
80  std::vector<edm::EDGetTokenT<QualityMaskCollection>> srcQuals;
81 
82  std::string priorityName_;
83 
84  float m_foundHitBonus;
85  float m_lostHitPenalty;
86  float m_shareFrac;
87  unsigned int m_minShareHits;
88  reco::TrackBase::TrackQuality m_minQuality;
89  bool m_allowFirstHitShare;
90  bool m_enableMerging;
91 
92  void produce(edm::StreamID, edm::Event& evt, const edm::EventSetup&) const override;
93 
94  bool areDuplicate(IHitV const& rh1, IHitV const& rh2) const;
95  };
96 
97 } // namespace
98 
100 
101 namespace {
102 
103  void TrackCollectionMerger::produce(edm::StreamID, edm::Event& evt, const edm::EventSetup& es) const {
104  TrackCollectionCloner::Producer producer(evt, collectionCloner);
105 
107  es.get<CkfComponentsRecord>().get(priorityName_, priorityH);
108  auto const& trackAlgoPriorityOrder = *priorityH;
109 
110  // load collections
111  auto collsSize = srcColls.size();
112  auto rSize = 0U;
113  declareDynArray(reco::TrackCollection const*, collsSize, trackColls);
114  for (auto i = 0U; i < collsSize; ++i) {
115  trackColls[i] = &srcColls[i].tracks(evt);
116  rSize += (*trackColls[i]).size();
117  }
118 
119  unsigned char qualMask = ~0;
120  const bool acceptAll = m_minQuality == reco::TrackBase::undefQuality;
121  if (!acceptAll)
122  qualMask = 1 << m_minQuality;
123 
124  // load tracks
125  initDynArray(unsigned int, collsSize, nGoods, 0);
126  declareDynArray(float, rSize, mvas);
127  declareDynArray(unsigned char, rSize, quals);
128  declareDynArray(unsigned int, rSize, tkInds);
129  auto k = 0U;
130  for (auto i = 0U; i < collsSize; ++i) {
131  auto const& tkColl = *trackColls[i];
132  auto size = tkColl.size();
134  evt.getByToken(srcMVAs[i], hmva);
135  assert((*hmva).size() == size);
137  evt.getByToken(srcQuals[i], hqual);
138  for (auto j = 0U; j < size; ++j) {
139  if (!(acceptAll || (qualMask & (*hqual)[j])))
140  continue;
141  mvas[k] = (*hmva)[j];
142  quals[k] = (*hqual)[j];
143  tkInds[k] = j;
144  ++k;
145  ++nGoods[i];
146  }
147  }
148 
149  auto ntotTk = k;
150  std::vector<bool> selected(ntotTk, true);
151 
154  declareDynArray(reco::TrackBase::AlgoMask, ntotTk, algoMask);
155 
156  // merge (if more than one collection...)
157 
158  auto merger = [&]() -> void {
159  // load momentum, hits and score
161  declareDynArray(float, ntotTk, score);
162  declareDynArray(IHitV, ntotTk, rh1);
163 
164  k = 0U;
165  for (auto i = 0U; i < collsSize; ++i) {
166  auto const& tkColl = *trackColls[i];
167  for (auto j = 0U; j < nGoods[i]; ++j) {
168  auto const& track = tkColl[tkInds[k]];
169  algo[k] = track.algo();
170  oriAlgo[k] = track.originalAlgo();
171  algoMask[k] = track.algoMask();
172  mom[k] = track.isLooper() ? reco::TrackBase::Vector() : track.momentum();
173  auto validHits = track.numberOfValidHits();
174  auto validPixelHits = track.hitPattern().numberOfValidPixelHits();
175  auto lostHits = track.numberOfLostHits();
176  score[k] = m_foundHitBonus * validPixelHits + m_foundHitBonus * validHits - m_lostHitPenalty * lostHits -
177  track.chi2();
178 
179  auto& rhv = rh1[k];
180  rhv.reserve(validHits);
181  auto compById = [](IHit const& h1, IHit const& h2) { return h1.first < h2.first; };
182  for (auto it = track.recHitsBegin(); it != track.recHitsEnd(); ++it) {
183  auto const& hit = *(*it);
184  auto id = hit.rawId();
185  if
186  LIKELY(hit.isValid()) {
187  rhv.emplace_back(id, &hit);
188  std::push_heap(rhv.begin(), rhv.end(), compById);
189  }
190  }
191  std::sort_heap(rhv.begin(), rhv.end(), compById);
192 
193  ++k;
194  }
195  }
196  assert(ntotTk == k);
197 
198  auto seti = [&](unsigned int ii, unsigned int jj) {
199  selected[jj] = false;
200  selected[ii] = true;
201  if (trackAlgoPriorityOrder.priority(oriAlgo[jj]) < trackAlgoPriorityOrder.priority(oriAlgo[ii]))
202  oriAlgo[ii] = oriAlgo[jj];
203  algoMask[ii] |= algoMask[jj];
204  quals[ii] |= (1 << reco::TrackBase::confirmed);
205  algoMask[jj] = algoMask[ii]; // in case we keep discarded
206  quals[jj] |= (1 << reco::TrackBase::discarded);
207  };
208 
209  auto iStart2 = 0U;
210  for (auto i = 0U; i < collsSize - 1; ++i) {
211  auto iStart1 = iStart2;
212  iStart2 = iStart1 + nGoods[i];
213  for (auto t1 = iStart1; t1 < iStart2; ++t1) {
214  if (!selected[t1])
215  continue;
216  auto score1 = score[t1];
217  for (auto t2 = iStart2; t2 < ntotTk; ++t2) {
218  if (!selected[t1])
219  break;
220  if (!selected[t2])
221  continue;
222  if (mom[t1].Dot(mom[t2]) < 0)
223  continue; // do not bother if in opposite hemespheres...
224  if (!areDuplicate(rh1[t1], rh1[t2]))
225  continue;
226  auto score2 = score[t2];
227 
228  // difference rather than ratio due to possible negative values for score
229  constexpr float almostSame = 0.01f;
230  if (score1 - score2 > almostSame) {
231  seti(t1, t2);
232  } else if (score2 - score1 > almostSame) {
233  seti(t2, t1);
234  } else { // take best
235  constexpr unsigned int qmask =
237  if ((quals[t1] & qmask) == (quals[t2] & qmask)) {
238  // take first
239  if (trackAlgoPriorityOrder.priority(algo[t1]) <= trackAlgoPriorityOrder.priority(algo[t2])) {
240  seti(t1, t2);
241  } else {
242  seti(t2, t1);
243  }
244  } else if ((quals[t1] & qmask) > (quals[t2] & qmask))
245  seti(t1, t2);
246  else
247  seti(t2, t1);
248  } // end ifs...
249 
250  } // end t2
251  } // end t1
252  } // end colls
253  }; // end merger;
254 
255  const bool doMerging = m_enableMerging && collsSize > 1;
256  if (doMerging)
257  merger();
258 
259  // products
260  auto pmvas = std::make_unique<MVACollection>();
261  auto pquals = std::make_unique<QualityMaskCollection>();
262 
263  // clone selected tracks...
264  auto nsel = 0U;
265  auto iStart2 = 0U;
266  auto isel = 0U;
267  for (auto i = 0U; i < collsSize; ++i) {
268  std::vector<unsigned int> selId;
269  std::vector<unsigned int> tid;
270  auto iStart1 = iStart2;
271  iStart2 = iStart1 + nGoods[i];
272  assert(producer.selTracks_->size() == isel);
273  for (auto t1 = iStart1; t1 < iStart2; ++t1) {
274  if (!selected[t1])
275  continue;
276  ++nsel;
277  tid.push_back(t1);
278  selId.push_back(tkInds[t1]);
279  pmvas->push_back(mvas[t1]);
280  pquals->push_back(quals[t1]);
281  }
282  producer(srcColls[i], selId);
283  assert(producer.selTracks_->size() == nsel);
284  assert(tid.size() == nsel - isel);
285  auto k = 0U;
286  for (; isel < nsel; ++isel) {
287  auto& otk = (*producer.selTracks_)[isel];
288  otk.setQualityMask((*pquals)[isel]); // needed also without merging
289  if (doMerging) {
290  otk.setOriginalAlgorithm(oriAlgo[tid[k]]);
291  otk.setAlgoMask(algoMask[tid[k++]]);
292  }
293  }
294  if (doMerging)
295  assert(tid.size() == k);
296  }
297 
298  assert(producer.selTracks_->size() == pmvas->size());
299 
300  // std::cout << "TrackCollectionMerger: sel tracks " << producer.selTracks_->size() << std::endl;
301 
302  evt.put(std::move(pmvas), "MVAValues");
303  evt.put(std::move(pquals), "QualityMasks");
304  }
305 
306  bool TrackCollectionMerger::areDuplicate(IHitV const& rh1, IHitV const& rh2) const {
307  auto nh1 = rh1.size();
308  auto nh2 = rh2.size();
309 
310  auto share = [](const TrackingRecHit* it, const TrackingRecHit* jt) -> bool {
311  return it->sharesInput(jt, TrackingRecHit::some);
312  };
313 
314  //loop over rechits
315  int noverlap = 0;
316  int firstoverlap = 0;
317  // check first hit (should use REAL first hit?)
318  if (m_allowFirstHitShare && rh1[0].first == rh2[0].first) {
319  if (share(rh1[0].second, rh2[0].second))
320  firstoverlap = 1;
321  }
322 
323  // exploit sorting
324  unsigned int jh = 0;
325  unsigned int ih = 0;
326  while (ih != nh1 && jh != nh2) {
327  auto const id1 = rh1[ih].first;
328  auto const id2 = rh2[jh].first;
329  if (id1 < id2)
330  ++ih;
331  else if (id2 < id1)
332  ++jh;
333  else {
334  if (share(rh1[ih].second, rh2[jh].second))
335  noverlap++;
336  ++jh;
337  ++ih;
338  } // equal ids
339  } //loop over ih & jh
340 
341  return noverlap >= int(m_minShareHits) &&
342  (noverlap - firstoverlap) > (std::min(nh1, nh2) - firstoverlap) * m_shareFrac;
343  }
344 
345 } // namespace
346 
349 
350 DEFINE_FWK_MODULE(TrackCollectionMerger);
ConfigurationDescriptions.h
trackAlgoPriorityOrder
constexpr std::array< unsigned int, reco::TrackBase::algoSize > trackAlgoPriorityOrder
Definition: trackAlgoPriorityOrder.h:157
edm::StreamID
Definition: StreamID.h:30
TrackExtra.h
RandomServiceHelper.t2
t2
Definition: RandomServiceHelper.py:257
electrons_cff.bool
bool
Definition: electrons_cff.py:372
mps_fire.i
i
Definition: mps_fire.py:355
edm::ParameterSetDescription::add
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Definition: ParameterSetDescription.h:95
TrackCollectionCloner::Producer
Definition: TrackCollectionCloner.h:52
reco::TrackBase::tight
Definition: TrackBase.h:153
reco::TrackBase::loose
Definition: TrackBase.h:152
ESHandle.h
min
T min(T a, T b)
Definition: MathUtil.h:58
reco::TrackBase::confirmed
Definition: TrackBase.h:155
reco::TrackBase::undefQuality
Definition: TrackBase.h:151
globals_cff.id1
id1
Definition: globals_cff.py:32
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
reco::TrackBase::TrackAlgorithm
TrackAlgorithm
track algorithm
Definition: TrackBase.h:89
TrackAlgoPriorityOrder.h
RawToDigi_Repacked_cff.producer
producer
Definition: RawToDigi_Repacked_cff.py:6
cms::cuda::assert
assert(be >=bs)
reco::TrackBase::TrackQuality
TrackQuality
track quality
Definition: TrackBase.h:150
edm::second
U second(std::pair< T, U > const &p)
Definition: ParameterSet.cc:215
edm::EDConsumerBase::consumesCollector
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
Definition: EDConsumerBase.cc:71
edm::global::EDProducerBase::produce
virtual void produce(StreamID, Event &, EventSetup const &) const =0
initDynArray
#define initDynArray(T, n, x, i)
Definition: DynArray.h:94
CkfComponentsRecord.h
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
edm::global::EDProducerBase::fillDescriptions
static void fillDescriptions(ConfigurationDescriptions &descriptions)
Definition: EDProducerBase.cc:255
edm::Handle
Definition: AssociativeIterator.h:50
dqmdumpme.first
first
Definition: dqmdumpme.py:55
TrackCollectionCloner
Definition: TrackCollectionCloner.h:25
cmsdt::algo
algo
Definition: constants.h:164
MakerMacros.h
TrackingRecHit.h
RandomServiceHelper.t1
t1
Definition: RandomServiceHelper.py:256
Track.h
edm::EventSetup::get
T get() const
Definition: EventSetup.h:73
TrackFwd.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
CkfComponentsRecord
Definition: CkfComponentsRecord.h:22
edm::ESHandle
Definition: DTSurvey.h:22
declareDynArray
#define declareDynArray(T, n, x)
Definition: DynArray.h:91
edm::Event::getByToken
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:528
dqmdumpme.k
k
Definition: dqmdumpme.py:60
ParameterSetDescription.h
TrackingRecHit::some
Definition: TrackingRecHit.h:59
TrajTrackAssociation.h
mitigatedMETSequence_cff.U
U
Definition: mitigatedMETSequence_cff.py:36
edm::global::EDProducer
Definition: EDProducer.h:32
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
TrackCollectionCloner::fill
static void fill(edm::ParameterSetDescription &desc)
Definition: TrackCollectionCloner.cc:3
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
funct::true
true
Definition: Factorize.h:173
edm::ParameterSet
Definition: ParameterSet.h:36
reco::TrackBase::discarded
Definition: TrackBase.h:159
Event.h
ModuleDef.h
createfilelist.int
int
Definition: createfilelist.py:10
edm::Event::put
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:132
edm::EventSetup
Definition: EventSetup.h:57
get
#define get
DynArray.h
TrackCollectionCloner.h
InputTag.h
TrackingRecHit
Definition: TrackingRecHit.h:21
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Trajectory.h
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
Frameworkfwd.h
electrons_cff.lostHits
lostHits
Definition: electrons_cff.py:407
LIKELY
#define LIKELY(x)
Definition: Likely.h:20
findQualityFiles.jj
string jj
Definition: findQualityFiles.py:188
reco::TrackBase::AlgoMask
std::bitset< algoSize > AlgoMask
algo mask
Definition: TrackBase.h:145
EventSetup.h
TrackingRecHit::sharesInput
virtual bool sharesInput(const TrackingRecHit *other, SharedInputType what) const
Definition: TrackingRecHit.cc:12
HLT_2018_cff.track
track
Definition: HLT_2018_cff.py:10352
ParameterSet.h
offlineSlimmedPrimaryVertices_cfi.score
score
Definition: offlineSlimmedPrimaryVertices_cfi.py:6
EDProducer.h
globals_cff.id2
id2
Definition: globals_cff.py:33
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
edm::Event
Definition: Event.h:73
cuy.ii
ii
Definition: cuy.py:590
edm::InputTag
Definition: InputTag.h:15
reco::TrackCollection
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
hit
Definition: SiStripHitEffFromCalibTree.cc:88
reco::TrackBase::Vector
math::XYZVector Vector
spatial vector
Definition: TrackBase.h:77
findQualityFiles.size
size
Write out results.
Definition: findQualityFiles.py:443
reco::TrackBase::highPurity
Definition: TrackBase.h:154