CMS 3D CMS Logo

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