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