CMS 3D CMS Logo

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