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