CMS 3D CMS Logo

TrackListMerger.cc
Go to the documentation of this file.
1 //
2 // Package: RecoTracker/FinalTrackSelectors
3 // Class: TrackListMerger
4 //
5 // Description: Hit Dumper
6 //
7 // Original Author: Steve Wagner, stevew@pizero.colorado.edu
8 // Created: Sat Jan 14 22:00:00 UTC 2006
9 //
10 //
11 
31 
33 public:
34  explicit TrackListMerger(const edm::ParameterSet& conf);
35 
36  ~TrackListMerger() override = default;
37 
38  void produce(edm::Event& e, const edm::EventSetup& c) override;
39 
40 private:
41  void returnEmptyCollections(edm::Event& e);
42 
43  using MVACollection = std::vector<float>;
44  using QualityMaskCollection = std::vector<unsigned char>;
45 
46  std::unique_ptr<reco::TrackCollection> outputTrks;
47  std::unique_ptr<reco::TrackExtraCollection> outputTrkExtras;
48  std::unique_ptr<TrackingRecHitCollection> outputTrkHits;
49  std::unique_ptr<std::vector<Trajectory>> outputTrajs;
50  std::unique_ptr<TrajTrackAssociationCollection> outputTTAss;
51  std::unique_ptr<TrajectorySeedCollection> outputSeeds;
52 
58 
61 
63 
64  struct TkEDGetTokenss {
74  edm::EDGetTokenT<std::vector<Trajectory>>&& traj_,
78  : tag(tag_), tk(tk_), traj(traj_), tass(tass_), tsel(tsel_), tmva(tmva_) {}
79  };
80  TkEDGetTokenss edTokens(const edm::InputTag& tag, const edm::InputTag& seltag, const edm::InputTag& mvatag) {
81  return TkEDGetTokenss(tag,
82  consumes<reco::TrackCollection>(tag),
83  consumes<std::vector<Trajectory>>(tag),
84  consumes<TrajTrackAssociationCollection>(tag),
85  consumes<edm::ValueMap<int>>(seltag),
86  consumes<edm::ValueMap<float>>(mvatag));
87  }
89  return TkEDGetTokenss(tag,
90  consumes<reco::TrackCollection>(tag),
91  consumes<std::vector<Trajectory>>(tag),
92  consumes<TrajTrackAssociationCollection>(tag),
94  consumes<edm::ValueMap<float>>(mvatag));
95  }
96  std::vector<TkEDGetTokenss> trackProducers_;
97 
99 
101  double minPT_;
102  unsigned int minFound_;
103  float epsilon_;
104  float shareFrac_;
107  std::vector<double> indivShareFrac_;
108 
109  std::vector<std::vector<int>> listsToMerge_;
110  std::vector<bool> promoteQuality_;
111  std::vector<int> hasSelector_;
112  bool copyMVA_;
113 
118 };
119 
120 #include <memory>
121 #include <string>
122 #include <iostream>
123 #include <cmath>
124 #include <vector>
125 
133 
135 
139 
140 //#include "DataFormats/TrackReco/src/classes.h"
141 
143 
144 #ifdef STAT_TSB
145 #include <x86intrin.h>
146 #endif
147 
148 namespace {
149 #ifdef STAT_TSB
150  inline volatile unsigned long long rdtsc() { return __rdtsc(); }
151 
152  struct StatCount {
153  float maxDP = 0.;
154  float maxDE = 0.;
155  unsigned long long st;
156  long long totBegin = 0;
157  long long totPre = 0;
158  long long totEnd = 0;
159  unsigned long long timeNo; // no-overlap
160  unsigned long long timeOv; // overlap
161  void begin(int tt) { totBegin += tt; }
162  void start() { st = rdtsc(); }
163  void noOverlap() { timeNo += (rdtsc() - st); }
164  void overlap() { timeOv += (rdtsc() - st); }
165  void pre(int tt) { totPre += tt; }
166  void end(int tt) { totEnd += tt; }
167  void de(float d) {
168  if (d > maxDE)
169  maxDE = d;
170  }
171  void dp(float d) {
172  if (d > maxDP)
173  maxDP = d;
174  }
175 
176  void print() const {
177  std::cout << "TrackListMerger stat\nBegin/Pre/End/maxDPhi/maxDEta/Overlap/NoOverlap " << totBegin << '/' << totPre
178  << '/' << totEnd << '/' << maxDP << '/' << maxDE << '/' << timeOv / 1000 << '/' << timeNo / 1000
179  << std::endl;
180  }
181  StatCount() {}
182  ~StatCount() { print(); }
183  };
184  StatCount statCount;
185 
186 #else
187  struct StatCount {
188  void begin(int) {}
189  void pre(int) {}
190  void end(int) {}
191  void start() {}
192  void noOverlap() {}
193  void overlap() {}
194  void de(float) {}
195  void dp(float) {}
196  };
197  CMS_THREAD_SAFE StatCount statCount;
198 #endif
199 
200 } // namespace
201 
202 namespace {
203  edm::ProductID clusterProductB(const TrackingRecHit* hit) {
204  return reinterpret_cast<const BaseTrackerRecHit*>(hit)->firstClusterRef().id();
205  }
206 } // namespace
207 
209  copyExtras_ = conf.getUntrackedParameter<bool>("copyExtras", true);
210  priorityName_ = conf.getParameter<std::string>("trackAlgoPriorityOrder");
211 
212  std::vector<edm::InputTag> trackProducerTags(conf.getParameter<std::vector<edm::InputTag>>("TrackProducers"));
213  //which of these do I need to turn into vectors?
214  maxNormalizedChisq_ = conf.getParameter<double>("MaxNormalizedChisq");
215  minPT_ = conf.getParameter<double>("MinPT");
216  minFound_ = (unsigned int)conf.getParameter<int>("MinFound");
217  epsilon_ = conf.getParameter<double>("Epsilon");
218  shareFrac_ = conf.getParameter<double>("ShareFrac");
219  allowFirstHitShare_ = conf.getParameter<bool>("allowFirstHitShare");
220  foundHitBonus_ = conf.getParameter<double>("FoundHitBonus");
221  lostHitPenalty_ = conf.getParameter<double>("LostHitPenalty");
222  indivShareFrac_ = conf.getParameter<std::vector<double>>("indivShareFrac");
223  std::string qualityStr = conf.getParameter<std::string>("newQuality");
224  priorityToken = esConsumes<TrackAlgoPriorityOrder, CkfComponentsRecord>(edm::ESInputTag("", priorityName_));
225 
226  if (!qualityStr.empty()) {
228  } else
230 
231  use_sharesInput_ = true;
232  if (epsilon_ > 0.0)
233  use_sharesInput_ = false;
234 
236 
237  for (unsigned int i = 0; i < setsToMerge.size(); i++) {
238  listsToMerge_.push_back(setsToMerge[i].getParameter<std::vector<int>>("tLists"));
239  promoteQuality_.push_back(setsToMerge[i].getParameter<bool>("pQual"));
240  }
241  hasSelector_ = conf.getParameter<std::vector<int>>("hasSelector");
242  copyMVA_ = conf.getParameter<bool>("copyMVA");
243 
244  std::vector<edm::InputTag> selectors(conf.getParameter<std::vector<edm::InputTag>>("selectedTrackQuals"));
245  std::vector<edm::InputTag> mvaStores;
246  if (conf.exists("mvaValueTags")) {
247  mvaStores = conf.getParameter<std::vector<edm::InputTag>>("mvaValueTags");
248  } else {
249  for (int i = 0; i < (int)selectors.size(); i++) {
250  edm::InputTag ntag(selectors[i].label(), "MVAVals");
251  mvaStores.push_back(ntag);
252  }
253  }
254  unsigned int numTrkColl = trackProducerTags.size();
255  if (numTrkColl != hasSelector_.size() || numTrkColl != selectors.size()) {
256  throw cms::Exception("Inconsistent size") << "need same number of track collections and selectors";
257  }
258  if (numTrkColl != hasSelector_.size() || numTrkColl != mvaStores.size()) {
259  throw cms::Exception("Inconsistent size") << "need same number of track collections and MVA stores";
260  }
261  for (unsigned int i = indivShareFrac_.size(); i < numTrkColl; i++) {
262  // edm::LogWarning("TrackListMerger") << "No indivShareFrac for " << trackProducersTags <<". Using default value of 1";
263  indivShareFrac_.push_back(1.0);
264  }
265 
266  trkQualMod_ = conf.getParameter<bool>("writeOnlyTrkQuals");
267  if (trkQualMod_) {
268  bool ok = true;
269  for (unsigned int i = 1; i < numTrkColl; i++) {
270  if (!(trackProducerTags[i] == trackProducerTags[0]))
271  ok = false;
272  }
273  if (!ok) {
274  throw cms::Exception("Bad input") << "to use writeOnlyTrkQuals=True all input InputTags must be the same";
275  }
276  produces<edm::ValueMap<int>>();
277  produces<QualityMaskCollection>("QualityMasks");
278  } else {
279  makeReKeyedSeeds_ = conf.getUntrackedParameter<bool>("makeReKeyedSeeds", false);
280  if (makeReKeyedSeeds_) {
281  copyExtras_ = true;
282  produces<TrajectorySeedCollection>();
283  }
284 
285  if (copyExtras_) {
286  produces<reco::TrackExtraCollection>();
287  produces<TrackingRecHitCollection>();
288  }
289 
290  // TrackCollection refers to TrackingRechit and TrackExtra
291  // collections, need to declare its production after them to work
292  // around a rare race condition in framework scheduling
293  produces<reco::TrackCollection>();
294 
295  produces<std::vector<Trajectory>>();
296  produces<TrajTrackAssociationCollection>();
297  }
298  produces<edm::ValueMap<float>>("MVAVals");
299  produces<MVACollection>("MVAValues");
300 
301  // Do all the consumes
302  trackProducers_.resize(numTrkColl);
303  for (unsigned int i = 0; i < numTrkColl; ++i) {
304  trackProducers_[i] = hasSelector_[i] > 0 ? edTokens(trackProducerTags[i], selectors[i], mvaStores[i])
305  : edTokens(trackProducerTags[i], mvaStores[i]);
306  }
307 }
308 
309 // Functions that gets called by framework every event
311  // extract tracker geometry
312  //
313  //edm::ESHandle<TrackerGeometry> theG;
314  //es.get<TrackerDigiGeometryRecord>().get(theG);
315 
316  // using namespace reco;
317 
319  auto const& trackAlgoPriorityOrder = *priorityH;
320 
321  // get Inputs
322  // if 1 input list doesn't exist, make an empty list, issue a warning, and continue
323  // this allows TrackListMerger to be used as a cleaner only if handed just one list
324  // if both input lists don't exist, will issue 2 warnings and generate an empty output collection
325  //
326  static const reco::TrackCollection s_empty;
327 
328  std::vector<const reco::TrackCollection*> trackColls;
329  std::vector<edm::Handle<reco::TrackCollection>> trackHandles(trackProducers_.size());
330  for (unsigned int i = 0; i < trackProducers_.size(); i++) {
331  trackColls.push_back(nullptr);
332  //edm::Handle<reco::TrackCollection> trackColl;
333  e.getByToken(trackProducers_[i].tk, trackHandles[i]);
334  if (trackHandles[i].isValid()) {
335  trackColls[i] = trackHandles[i].product();
336  } else {
337  edm::LogWarning("TrackListMerger") << "TrackCollection " << trackProducers_[i].tag << " not found";
338  trackColls[i] = &s_empty;
339  }
340  }
341 
342  if (trackColls.empty()) {
343  // output empty collections and early return
344  this->returnEmptyCollections(e);
345  return;
346  }
347 
348  unsigned int collsSize = trackColls.size();
349  unsigned int rSize = 0;
350  unsigned int trackCollSizes[collsSize];
351  unsigned int trackCollFirsts[collsSize];
352  for (unsigned int i = 0; i != collsSize; i++) {
353  trackCollSizes[i] = trackColls[i]->size();
354  trackCollFirsts[i] = rSize;
355  rSize += trackCollSizes[i];
356  }
357 
358  if (rSize == 0) {
359  // output empty collections and early return
360  this->returnEmptyCollections(e);
361  return;
362  }
363 
364  statCount.begin(rSize);
365 
366  //
367  // quality cuts first
368  //
369  int i = -1;
370 
371  int selected[rSize];
372  int indexG[rSize];
373  bool trkUpdated[rSize];
374  int trackCollNum[rSize];
375  int trackQuals[rSize];
376  float trackMVAs[rSize];
377  reco::TrackBase::TrackAlgorithm oriAlgo[rSize];
378  std::vector<reco::TrackBase::AlgoMask> algoMask(rSize);
379  for (unsigned int j = 0; j < rSize; j++) {
380  indexG[j] = -1;
381  selected[j] = 1;
382  trkUpdated[j] = false;
383  trackCollNum[j] = 0;
384  trackQuals[j] = 0;
385  trackMVAs[j] = -998.0;
387  }
388 
389  int ngood = 0;
390  for (unsigned int j = 0; j != collsSize; j++) {
391  const reco::TrackCollection* tC1 = trackColls[j];
392 
393  edm::Handle<edm::ValueMap<int>> trackSelColl;
394  edm::Handle<edm::ValueMap<float>> trackMVAStore;
395  if (copyMVA_)
396  e.getByToken(trackProducers_[j].tmva, trackMVAStore);
397  if (hasSelector_[j] > 0) {
398  e.getByToken(trackProducers_[j].tsel, trackSelColl);
399  }
400 
401  if (!tC1->empty()) {
402  unsigned int iC = 0;
403  for (reco::TrackCollection::const_iterator track = tC1->begin(); track != tC1->end(); track++) {
404  i++;
405  trackCollNum[i] = j;
406  trackQuals[i] = track->qualityMask();
407  oriAlgo[i] = track->originalAlgo();
408  algoMask[i] = track->algoMask();
409 
410  reco::TrackRef trkRef = reco::TrackRef(trackHandles[j], iC);
411  if (copyMVA_)
412  if ((*trackMVAStore).contains(trkRef.id()))
413  trackMVAs[i] = (*trackMVAStore)[trkRef];
414  if (hasSelector_[j] > 0) {
415  int qual = (*trackSelColl)[trkRef];
416  if (qual < 0) {
417  selected[i] = 0;
418  iC++;
419  continue;
420  } else {
421  trackQuals[i] = qual;
422  }
423  }
424  iC++;
425  selected[i] = trackQuals[i] + 10; //10 is magic number used throughout...
426  if ((short unsigned)track->ndof() < 1) {
427  selected[i] = 0;
428  continue;
429  }
430  if (track->normalizedChi2() > maxNormalizedChisq_) {
431  selected[i] = 0;
432  continue;
433  }
434  if (track->found() < minFound_) {
435  selected[i] = 0;
436  continue;
437  }
438  if (track->pt() < minPT_) {
439  selected[i] = 0;
440  continue;
441  }
442  // good!
443  indexG[i] = ngood++;
444  //if ( beVerb) std::cout << "inverb " << track->pt() << " " << selected[i] << std::endl;
445  } //end loop over tracks
446  } //end more than 0 track
447  } // loop over trackcolls
448 
449  statCount.pre(ngood);
450 
451  //cache the id and rechits of valid hits
452  typedef std::pair<unsigned int, const TrackingRecHit*> IHit;
453  std::vector<std::vector<IHit>> rh1(ngood); // "not an array" of vectors!
454  //const TrackingRecHit* fh1[ngood]; // first hit...
456  float score[std::max(1, ngood)];
457 
458  for (unsigned int j = 0; j < rSize; j++) {
459  if (selected[j] == 0)
460  continue;
461  int i = indexG[j];
462  assert(i >= 0);
463  unsigned int collNum = trackCollNum[j];
464  unsigned int trackNum = j - trackCollFirsts[collNum];
465  const reco::Track* track = &((trackColls[collNum])->at(trackNum));
466 
467  algo[i] = track->algo();
468  int validHits = track->numberOfValidHits();
469  int validPixelHits = track->hitPattern().numberOfValidPixelHits();
470  int lostHits = track->numberOfLostHits();
471  score[i] =
472  foundHitBonus_ * validPixelHits + foundHitBonus_ * validHits - lostHitPenalty_ * lostHits - track->chi2();
473 
474  rh1[i].reserve(validHits);
475  auto compById = [](IHit const& h1, IHit const& h2) { return h1.first < h2.first; };
476  for (trackingRecHit_iterator it = track->recHitsBegin(); it != track->recHitsEnd(); ++it) {
477  const TrackingRecHit* hit = (*it);
478  unsigned int id = hit->rawId();
479  if (hit->geographicalId().subdetId() > 2)
480  id &= (~3); // mask mono/stereo in strips...
481  if LIKELY (hit->isValid()) {
482  rh1[i].emplace_back(id, hit);
483  std::push_heap(rh1[i].begin(), rh1[i].end(), compById);
484  }
485  }
486  std::sort_heap(rh1[i].begin(), rh1[i].end(), compById);
487  }
488 
489  //DL here
490  if LIKELY (ngood > 1 && collsSize > 1)
491  for (unsigned int ltm = 0; ltm < listsToMerge_.size(); ltm++) {
492  int saveSelected[rSize];
493  bool notActive[collsSize];
494  for (unsigned int cn = 0; cn != collsSize; ++cn)
495  notActive[cn] = find(listsToMerge_[ltm].begin(), listsToMerge_[ltm].end(), cn) == listsToMerge_[ltm].end();
496 
497  for (unsigned int i = 0; i < rSize; i++)
498  saveSelected[i] = selected[i];
499 
500  //DL protect against 0 tracks?
501  for (unsigned int i = 0; i < rSize - 1; i++) {
502  if (selected[i] == 0)
503  continue;
504  unsigned int collNum = trackCollNum[i];
505 
506  //check that this track is in one of the lists for this iteration
507  if (notActive[collNum])
508  continue;
509 
510  int k1 = indexG[i];
511  unsigned int nh1 = rh1[k1].size();
512  int qualityMaskT1 = trackQuals[i];
513 
514  int nhit1 = nh1; // validHits[k1];
515  float score1 = score[k1];
516 
517  // start at next collection
518  for (unsigned int j = i + 1; j < rSize; j++) {
519  if (selected[j] == 0)
520  continue;
521  unsigned int collNum2 = trackCollNum[j];
522  if ((collNum == collNum2) && indivShareFrac_[collNum] > 0.99)
523  continue;
524  //check that this track is in one of the lists for this iteration
525  if (notActive[collNum2])
526  continue;
527 
528  int k2 = indexG[j];
529 
530  int newQualityMask = -9; //avoid resetting quality mask if not desired 10+ -9 =1
531  if (promoteQuality_[ltm]) {
532  int maskT1 = saveSelected[i] > 1 ? saveSelected[i] - 10 : qualityMaskT1;
533  int maskT2 = saveSelected[j] > 1 ? saveSelected[j] - 10 : trackQuals[j];
534  newQualityMask = (maskT1 | maskT2); // take OR of trackQuality
535  }
536  unsigned int nh2 = rh1[k2].size();
537  int nhit2 = nh2;
538 
539  auto share = use_sharesInput_ ? [](const TrackingRecHit* it, const TrackingRecHit* jt, float) -> bool {
540  return it->sharesInput(jt, TrackingRecHit::some);
541  }
542  : [](const TrackingRecHit* it, const TrackingRecHit* jt, float eps) -> bool {
543  float delta = std::abs(it->localPosition().x() - jt->localPosition().x());
544  return (it->geographicalId() == jt->geographicalId()) && (delta < eps);
545  };
546 
547  statCount.start();
548 
549  //loop over rechits
550  int noverlap = 0;
551  int firstoverlap = 0;
552  // check first hit (should use REAL first hit?)
553  if UNLIKELY (allowFirstHitShare_ && rh1[k1][0].first == rh1[k2][0].first) {
554  const TrackingRecHit* it = rh1[k1][0].second;
555  const TrackingRecHit* jt = rh1[k2][0].second;
556  if (share(it, jt, epsilon_))
557  firstoverlap = 1;
558  }
559 
560  // exploit sorting
561  unsigned int jh = 0;
562  unsigned int ih = 0;
563  while (ih != nh1 && jh != nh2) {
564  // break if not enough to go...
565  // if ( nprecut-noverlap+firstoverlap > int(nh1-ih)) break;
566  // if ( nprecut-noverlap+firstoverlap > int(nh2-jh)) break;
567  auto const id1 = rh1[k1][ih].first;
568  auto const id2 = rh1[k2][jh].first;
569  if (id1 < id2)
570  ++ih;
571  else if (id2 < id1)
572  ++jh;
573  else {
574  // in case of split-hit do full conbinatorics
575  auto li = ih;
576  while ((++li) != nh1 && id1 == rh1[k1][li].first) {
577  }
578  auto lj = jh;
579  while ((++lj) != nh2 && id2 == rh1[k2][lj].first) {
580  }
581  for (auto ii = ih; ii != li; ++ii)
582  for (auto jj = jh; jj != lj; ++jj) {
583  const TrackingRecHit* it = rh1[k1][ii].second;
584  const TrackingRecHit* jt = rh1[k2][jj].second;
585  if (share(it, jt, epsilon_))
586  noverlap++;
587  }
588  jh = lj;
589  ih = li;
590  } // equal ids
591 
592  } //loop over ih & jh
593 
594  bool dupfound =
595  (collNum != collNum2)
596  ? (noverlap - firstoverlap) > (std::min(nhit1, nhit2) - firstoverlap) * shareFrac_
597  : (noverlap - firstoverlap) > (std::min(nhit1, nhit2) - firstoverlap) * indivShareFrac_[collNum];
598 
599  auto seti = [&](unsigned int ii, unsigned int jj) {
600  selected[jj] = 0;
601  selected[ii] = 10 + newQualityMask; // add 10 to avoid the case where mask = 1
602  trkUpdated[ii] = true;
603  if (trackAlgoPriorityOrder.priority(oriAlgo[jj]) < trackAlgoPriorityOrder.priority(oriAlgo[ii]))
604  oriAlgo[ii] = oriAlgo[jj];
605  algoMask[ii] |= algoMask[jj];
606  algoMask[jj] = algoMask[ii]; // in case we keep discarded
607  };
608 
609  if (dupfound) {
610  float score2 = score[k2];
611  constexpr float almostSame =
612  0.01f; // difference rather than ratio due to possible negative values for score
613  if (score1 - score2 > almostSame) {
614  seti(i, j);
615  } else if (score2 - score1 > almostSame) {
616  seti(j, i);
617  } else {
618  // If tracks from both iterations are virtually identical, choose the one with the best quality or with lower algo
619  if ((trackQuals[j] &
621  (trackQuals[i] &
623  //same quality, pick earlier algo
624  if (trackAlgoPriorityOrder.priority(algo[k1]) <= trackAlgoPriorityOrder.priority(algo[k2])) {
625  seti(i, j);
626  } else {
627  seti(j, i);
628  }
629  } else if ((trackQuals[j] & (1 << reco::TrackBase::loose | 1 << reco::TrackBase::tight |
631  (trackQuals[i] & (1 << reco::TrackBase::loose | 1 << reco::TrackBase::tight |
633  seti(i, j);
634  } else {
635  seti(j, i);
636  }
637  } //end fi < fj
638  statCount.overlap();
639  /*
640  if (at0[k1]&&at0[k2]) {
641  statCount.dp(dphi);
642  if (dz<1.f) statCount.de(deta);
643  }
644  */
645  } //end got a duplicate
646  else {
647  statCount.noOverlap();
648  }
649  //stop if the ith track is now unselected
650  if (selected[i] == 0)
651  break;
652  } //end track2 loop
653  } //end track loop
654  } //end loop over track list sets
655 
656  auto vmMVA = std::make_unique<edm::ValueMap<float>>();
657  edm::ValueMap<float>::Filler fillerMVA(*vmMVA);
658 
659  // special case - if just doing the trkquals
660  if (trkQualMod_) {
661  unsigned int tSize = trackColls[0]->size();
662  auto vm = std::make_unique<edm::ValueMap<int>>();
664 
665  std::vector<int> finalQuals(tSize, -1); //default is unselected
666  for (unsigned int i = 0; i < rSize; i++) {
667  unsigned int tNum = i % tSize;
668 
669  if (selected[i] > 1) {
670  finalQuals[tNum] = selected[i] - 10;
671  if (trkUpdated[i])
672  finalQuals[tNum] = (finalQuals[tNum] | (1 << qualityToSet_));
673  }
674  if (selected[i] == 1)
675  finalQuals[tNum] = trackQuals[i];
676  }
677 
678  filler.insert(trackHandles[0], finalQuals.begin(), finalQuals.end());
679  filler.fill();
680 
681  e.put(std::move(vm));
682  for (auto& q : finalQuals)
683  q = std::max(q, 0);
684  auto quals = std::make_unique<QualityMaskCollection>(finalQuals.begin(), finalQuals.end());
685  e.put(std::move(quals), "QualityMasks");
686 
687  std::vector<float> mvaVec(tSize, -99);
688 
689  for (unsigned int i = 0; i < rSize; i++) {
690  unsigned int tNum = i % tSize;
691  mvaVec[tNum] = trackMVAs[tNum];
692  }
693 
694  fillerMVA.insert(trackHandles[0], mvaVec.begin(), mvaVec.end());
695  fillerMVA.fill();
696  if (copyMVA_) {
697  e.put(std::move(vmMVA), "MVAVals");
698  auto mvas = std::make_unique<MVACollection>(mvaVec.begin(), mvaVec.end());
699  e.put(std::move(mvas), "MVAValues");
700  }
701  return;
702  }
703 
704  //
705  // output selected tracks - if any
706  //
707 
708  std::vector<reco::TrackRef> trackRefs(rSize);
709  std::vector<edm::RefToBase<TrajectorySeed>> seedsRefs(rSize);
710 
711  unsigned int nToWrite = 0;
712  for (unsigned int i = 0; i < rSize; i++)
713  if (selected[i] != 0)
714  nToWrite++;
715 
716  std::vector<float> mvaVec;
717 
718  outputTrks = std::make_unique<reco::TrackCollection>();
719  outputTrks->reserve(nToWrite);
720  refTrks = e.getRefBeforePut<reco::TrackCollection>();
721 
722  if (copyExtras_) {
723  outputTrkExtras = std::make_unique<reco::TrackExtraCollection>();
724  outputTrkExtras->reserve(nToWrite);
725  refTrkExtras = e.getRefBeforePut<reco::TrackExtraCollection>();
726  outputTrkHits = std::make_unique<TrackingRecHitCollection>();
727  outputTrkHits->reserve(nToWrite * 25);
728  refTrkHits = e.getRefBeforePut<TrackingRecHitCollection>();
729  if (makeReKeyedSeeds_) {
730  outputSeeds = std::make_unique<TrajectorySeedCollection>();
731  outputSeeds->reserve(nToWrite);
732  refTrajSeeds = e.getRefBeforePut<TrajectorySeedCollection>();
733  }
734  }
735 
736  outputTrajs = std::make_unique<std::vector<Trajectory>>();
737  outputTrajs->reserve(rSize);
738 
739  for (unsigned int i = 0; i < rSize; i++) {
740  if (selected[i] == 0) {
741  trackRefs[i] = reco::TrackRef();
742  continue;
743  }
744 
745  unsigned int collNum = trackCollNum[i];
746  unsigned int trackNum = i - trackCollFirsts[collNum];
747  const reco::Track* track = &((trackColls[collNum])->at(trackNum));
748  outputTrks->push_back(reco::Track(*track));
749  mvaVec.push_back(trackMVAs[i]);
750  if (selected[i] > 1) {
751  outputTrks->back().setQualityMask(selected[i] - 10);
752  if (trkUpdated[i])
753  outputTrks->back().setQuality(qualityToSet_);
754  }
755  //might duplicate things, but doesnt hurt
756  if (selected[i] == 1)
757  outputTrks->back().setQualityMask(trackQuals[i]);
758  outputTrks->back().setOriginalAlgorithm(oriAlgo[i]);
759  outputTrks->back().setAlgoMask(algoMask[i]);
760 
761  // if ( beVerb ) std::cout << "selected " << outputTrks->back().pt() << " " << outputTrks->back().qualityMask() << " " << selected[i] << std::endl;
762 
763  //fill the TrackCollection
764  if (copyExtras_) {
765  edm::RefToBase<TrajectorySeed> origSeedRef = track->seedRef();
766  //creating a seed with rekeyed clusters if required
767  if (makeReKeyedSeeds_) {
768  bool doRekeyOnThisSeed = false;
769 
770  edm::InputTag clusterRemovalInfos("");
771  //grab on of the hits of the seed
772  if (origSeedRef->nHits() != 0) {
773  TrackingRecHit const& hit = *origSeedRef->recHits().begin();
774  if (hit.isValid()) {
775  edm::ProductID pID = clusterProductB(&hit);
776  // the cluster collection either produced a removalInfo or mot
777  //get the clusterremoval info from the provenance: will rekey if this is found
779  edm::StableProvenance const& prov = e.getStableProvenance(pID);
780  clusterRemovalInfos = edm::InputTag(prov.moduleLabel(), prov.productInstanceName(), prov.processName());
781  doRekeyOnThisSeed = e.getByLabel(clusterRemovalInfos, CRIh);
782  } //valid hit
783  } //nhit!=0
784 
785  if (doRekeyOnThisSeed && !(clusterRemovalInfos == edm::InputTag(""))) {
786  ClusterRemovalRefSetter refSetter(e, clusterRemovalInfos);
787  TrajectorySeed::RecHitContainer newRecHitContainer;
788  newRecHitContainer.reserve(origSeedRef->nHits());
789  for (auto const& recHit : origSeedRef->recHits()) {
790  newRecHitContainer.push_back(recHit);
791  refSetter.reKey(&newRecHitContainer.back());
792  }
793  outputSeeds->push_back(
794  TrajectorySeed(origSeedRef->startingState(), newRecHitContainer, origSeedRef->direction()));
795  }
796  //doRekeyOnThisSeed=true
797  else {
798  //just copy the one we had before
799  outputSeeds->push_back(TrajectorySeed(*origSeedRef));
800  }
802  origSeedRef = edm::RefToBase<TrajectorySeed>(pureRef);
803  } //creating a new seed and rekeying it rechit clusters.
804 
805  // Fill TrackExtra collection
806  outputTrkExtras->push_back(reco::TrackExtra(track->outerPosition(),
807  track->outerMomentum(),
808  track->outerOk(),
809  track->innerPosition(),
810  track->innerMomentum(),
811  track->innerOk(),
812  track->outerStateCovariance(),
813  track->outerDetId(),
814  track->innerStateCovariance(),
815  track->innerDetId(),
816  track->seedDirection(),
817  origSeedRef));
818  seedsRefs[i] = origSeedRef;
819  outputTrks->back().setExtra(reco::TrackExtraRef(refTrkExtras, outputTrkExtras->size() - 1));
820  reco::TrackExtra& tx = outputTrkExtras->back();
821  tx.setResiduals(track->residuals());
822 
823  // fill TrackingRecHits
824  unsigned nh1 = track->recHitsSize();
825  tx.setHits(refTrkHits, outputTrkHits->size(), nh1);
826  tx.setTrajParams(track->extra()->trajParams(), track->extra()->chi2sX5());
827  assert(tx.trajParams().size() == tx.recHitsSize());
828  for (auto hh = track->recHitsBegin(), eh = track->recHitsEnd(); hh != eh; ++hh) {
829  outputTrkHits->push_back((*hh)->clone());
830  }
831  }
832  trackRefs[i] = reco::TrackRef(refTrks, outputTrks->size() - 1);
833 
834  } //end faux loop over tracks
835 
836  //Fill the trajectories, etc. for 1st collection
837  refTrajs = e.getRefBeforePut<std::vector<Trajectory>>();
838 
839  outputTTAss = std::make_unique<TrajTrackAssociationCollection>(refTrajs, refTrks);
840 
841  for (unsigned int ti = 0; ti < trackColls.size(); ti++) {
844  e.getByToken(trackProducers_[ti].traj, hTraj1);
845  e.getByToken(trackProducers_[ti].tass, hTTAss1);
846 
847  if (hTraj1.failedToGet() || hTTAss1.failedToGet())
848  continue;
849 
850  for (size_t i = 0, n = hTraj1->size(); i < n; ++i) {
851  edm::Ref<std::vector<Trajectory>> trajRef(hTraj1, i);
853  if (match != hTTAss1->end()) {
854  const edm::Ref<reco::TrackCollection>& trkRef = match->val;
855  uint32_t oldKey = trackCollFirsts[ti] + static_cast<uint32_t>(trkRef.key());
856  if (trackRefs[oldKey].isNonnull()) {
857  outputTrajs->push_back(*trajRef);
858  //if making extras and the seeds at the same time, change the seed ref on the trajectory
860  outputTrajs->back().setSeedRef(seedsRefs[oldKey]);
861  outputTTAss->insert(edm::Ref<std::vector<Trajectory>>(refTrajs, outputTrajs->size() - 1), trackRefs[oldKey]);
862  }
863  }
864  }
865  }
866 
867  statCount.end(outputTrks->size());
868 
869  edm::ProductID nPID = refTrks.id();
870  edm::TestHandle<reco::TrackCollection> outHandle(outputTrks.get(), nPID);
871  fillerMVA.insert(outHandle, mvaVec.begin(), mvaVec.end());
872  fillerMVA.fill();
873 
874  e.put(std::move(outputTrks));
875  if (copyMVA_) {
876  e.put(std::move(vmMVA), "MVAVals");
877  auto mvas = std::make_unique<MVACollection>(mvaVec.begin(), mvaVec.end());
878  e.put(std::move(mvas), "MVAValues");
879  }
880  if (copyExtras_) {
882  e.put(std::move(outputTrkHits));
883  if (makeReKeyedSeeds_)
884  e.put(std::move(outputSeeds));
885  }
886  e.put(std::move(outputTrajs));
887  e.put(std::move(outputTTAss));
888  return;
889 
890 } //end produce
891 
893  if (trkQualMod_) {
894  auto vm = std::make_unique<edm::ValueMap<int>>();
895  e.put(std::move(vm));
896  auto quals = std::make_unique<QualityMaskCollection>();
897  e.put(std::move(quals), "QualityMasks");
898  } else {
899  auto outputTrks = std::make_unique<reco::TrackCollection>();
900  e.put(std::move(outputTrks));
901 
902  if (makeReKeyedSeeds_) {
903  auto outputSeeds = std::make_unique<TrajectorySeedCollection>();
904  e.put(std::move(outputSeeds));
905  }
906 
907  if (copyExtras_) {
908  auto outputTrkExtras = std::make_unique<reco::TrackExtraCollection>();
909  auto outputTrkHits = std::make_unique<TrackingRecHitCollection>();
911  e.put(std::move(outputTrkHits));
912  }
913 
914  auto outputTrajs = std::make_unique<std::vector<Trajectory>>();
915  outputTTAss = std::make_unique<TrajTrackAssociationCollection>();
916  e.put(std::move(outputTrajs));
917  e.put(std::move(outputTTAss));
918  }
919  auto vmMVA = std::make_unique<edm::ValueMap<float>>();
920  e.put(std::move(vmMVA), "MVAVals");
921  auto mvas = std::make_unique<MVACollection>();
922  e.put(std::move(mvas), "MVAValues");
923  return;
924 }
925 
928 
TrackListMerger(const edm::ParameterSet &conf)
std::vector< unsigned char > QualityMaskCollection
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
reference back()
Definition: OwnVector.h:431
TrajParams const & trajParams() const
std::vector< bool > promoteQuality_
std::unique_ptr< TrajTrackAssociationCollection > outputTTAss
ProductID id() const
Accessor for product ID.
Definition: Ref.h:244
const bool isValid(const Frame &aFrame, const FrameQuality &aQuality, const uint16_t aExpectedPos)
T begin() const
Definition: Range.h:15
reco::TrackRefProd refTrks
RecHitRange recHits() const
void returnEmptyCollections(edm::Event &e)
TkEDGetTokenss edTokens(const edm::InputTag &tag, const edm::InputTag &mvatag)
TrackQuality
track quality
Definition: TrackBase.h:150
void setHits(TrackingRecHitRefProd const &prod, unsigned firstH, unsigned int nH)
std::string const & productInstanceName() const
std::vector< ParameterSet > VParameterSet
Definition: ParameterSet.h:35
unsigned int recHitsSize() const
number of RecHits
bool exists(std::string const &parameterName) const
checks if a parameter exists
reco::TrackBase::TrackQuality qualityToSet_
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
#define LIKELY(x)
Definition: Likely.h:20
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
assert(be >=bs)
std::vector< float > MVACollection
std::unique_ptr< std::vector< Trajectory > > outputTrajs
key_type key() const
Accessor for product key.
Definition: Ref.h:250
const_iterator find(const key_type &k) const
find element with specified reference key
TrackAlgorithm
track algorithm
Definition: TrackBase.h:89
TrackingRecHitRefProd refTrkHits
void produce(edm::Event &e, const edm::EventSetup &c) override
std::vector< double > indivShareFrac_
const_iterator end() const
last iterator over the map (read only)
T getUntrackedParameter(std::string const &, T const &) const
PropagationDirection direction() const
bool failedToGet() const
Definition: HandleBase.h:72
char const * label
void push_back(D *&d)
Definition: OwnVector.h:326
std::unique_ptr< TrackingRecHitCollection > outputTrkHits
Definition: TTTypes.h:54
edm::EDGetTokenT< edm::ValueMap< int > > tsel
std::vector< TrajectorySeed > TrajectorySeedCollection
std::string const & moduleLabel() const
unsigned int nHits() const
edm::RefProd< TrajectorySeedCollection > refTrajSeeds
edm::EDGetTokenT< reco::TrackCollection > tk
void print(TMatrixD &m, const char *label=nullptr, bool mathematicaFormat=false)
Definition: Utilities.cc:47
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::string const & processName() const
TkEDGetTokenss edTokens(const edm::InputTag &tag, const edm::InputTag &seltag, const edm::InputTag &mvatag)
#define CMS_THREAD_SAFE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::ESGetToken< TrackAlgoPriorityOrder, CkfComponentsRecord > priorityToken
edm::EDGetTokenT< TrajTrackAssociationCollection > tass
edm::RefProd< std::vector< Trajectory > > refTrajs
ALPAKA_FN_ACC ALPAKA_FN_INLINE void uint32_t const uint32_t CACellT< TrackerTraits > uint32_t CellNeighborsVector< TrackerTraits > CellTracksVector< TrackerTraits > HitsConstView< TrackerTraits > hh
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
d
Definition: ztail.py:151
PTrajectoryStateOnDet const & startingState() const
std::unique_ptr< reco::TrackExtraCollection > outputTrkExtras
std::vector< int > hasSelector_
ii
Definition: cuy.py:589
ProductID id() const
Accessor for product ID.
Definition: RefProd.h:124
#define dso_hidden
Definition: Visibility.h:12
std::vector< TrackExtra > TrackExtraCollection
collection of TrackExtra objects
Definition: TrackExtraFwd.h:10
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:126
std::unique_ptr< TrajectorySeedCollection > outputSeeds
edm::EDGetTokenT< edm::ValueMap< float > > tmva
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:20
constexpr std::array< unsigned int, reco::TrackBase::algoSize > trackAlgoPriorityOrder
void reKey(TrackingRecHit *hit) const
unsigned int minFound_
std::vector< TkEDGetTokenss > trackProducers_
TkEDGetTokenss(const edm::InputTag &tag_, edm::EDGetTokenT< reco::TrackCollection > &&tk_, edm::EDGetTokenT< std::vector< Trajectory >> &&traj_, edm::EDGetTokenT< TrajTrackAssociationCollection > &&tass_, edm::EDGetTokenT< edm::ValueMap< int >> &&tsel_, edm::EDGetTokenT< edm::ValueMap< float >> &&tmva_)
std::vector< std::vector< int > > listsToMerge_
std::unique_ptr< reco::TrackCollection > outputTrks
edm::EDGetTokenT< std::vector< Trajectory > > traj
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
#define UNLIKELY(x)
Definition: Likely.h:21
Log< level::Warning, false > LogWarning
std::string priorityName_
def move(src, dest)
Definition: eostools.py:511
void reserve(size_t)
Definition: OwnVector.h:320
void setTrajParams(TrajParams tmps, Chi2sFive chi2s)
reco::TrackExtraRefProd refTrkExtras
void setResiduals(const TrackResiduals &r)
set the residuals
Definition: TrackExtra.h:132