CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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;
37 
38  void produce(edm::Event& e, const edm::EventSetup& c) override;
39 
40 private:
41  using MVACollection = std::vector<float>;
42  using QualityMaskCollection = std::vector<unsigned char>;
43 
44  std::unique_ptr<reco::TrackCollection> outputTrks;
45  std::unique_ptr<reco::TrackExtraCollection> outputTrkExtras;
46  std::unique_ptr<TrackingRecHitCollection> outputTrkHits;
47  std::unique_ptr<std::vector<Trajectory>> outputTrajs;
48  std::unique_ptr<TrajTrackAssociationCollection> outputTTAss;
49  std::unique_ptr<TrajectorySeedCollection> outputSeeds;
50 
56 
59 
61 
62  struct TkEDGetTokenss {
72  edm::EDGetTokenT<std::vector<Trajectory>>&& traj_,
76  : tag(tag_), tk(tk_), traj(traj_), tass(tass_), tsel(tsel_), tmva(tmva_) {}
77  };
78  TkEDGetTokenss edTokens(const edm::InputTag& tag, const edm::InputTag& seltag, const edm::InputTag& mvatag) {
79  return TkEDGetTokenss(tag,
80  consumes<reco::TrackCollection>(tag),
81  consumes<std::vector<Trajectory>>(tag),
82  consumes<TrajTrackAssociationCollection>(tag),
83  consumes<edm::ValueMap<int>>(seltag),
84  consumes<edm::ValueMap<float>>(mvatag));
85  }
87  return TkEDGetTokenss(tag,
88  consumes<reco::TrackCollection>(tag),
89  consumes<std::vector<Trajectory>>(tag),
90  consumes<TrajTrackAssociationCollection>(tag),
92  consumes<edm::ValueMap<float>>(mvatag));
93  }
94  std::vector<TkEDGetTokenss> trackProducers_;
95 
97 
99  double minPT_;
100  unsigned int minFound_;
101  float epsilon_;
102  float shareFrac_;
105  std::vector<double> indivShareFrac_;
106 
107  std::vector<std::vector<int>> listsToMerge_;
108  std::vector<bool> promoteQuality_;
109  std::vector<int> hasSelector_;
110  bool copyMVA_;
111 
116 };
117 
118 #include <memory>
119 #include <string>
120 #include <iostream>
121 #include <cmath>
122 #include <vector>
123 
131 
133 
137 
138 //#include "DataFormats/TrackReco/src/classes.h"
139 
141 
142 #ifdef STAT_TSB
143 #include <x86intrin.h>
144 #endif
145 
146 namespace {
147 #ifdef STAT_TSB
148  inline volatile unsigned long long rdtsc() { return __rdtsc(); }
149 
150  struct StatCount {
151  float maxDP = 0.;
152  float maxDE = 0.;
153  unsigned long long st;
154  long long totBegin = 0;
155  long long totPre = 0;
156  long long totEnd = 0;
157  unsigned long long timeNo; // no-overlap
158  unsigned long long timeOv; // overlap
159  void begin(int tt) { totBegin += tt; }
160  void start() { st = rdtsc(); }
161  void noOverlap() { timeNo += (rdtsc() - st); }
162  void overlap() { timeOv += (rdtsc() - st); }
163  void pre(int tt) { totPre += tt; }
164  void end(int tt) { totEnd += tt; }
165  void de(float d) {
166  if (d > maxDE)
167  maxDE = d;
168  }
169  void dp(float d) {
170  if (d > maxDP)
171  maxDP = d;
172  }
173 
174  void print() const {
175  std::cout << "TrackListMerger stat\nBegin/Pre/End/maxDPhi/maxDEta/Overlap/NoOverlap " << totBegin << '/' << totPre
176  << '/' << totEnd << '/' << maxDP << '/' << maxDE << '/' << timeOv / 1000 << '/' << timeNo / 1000
177  << std::endl;
178  }
179  StatCount() {}
180  ~StatCount() { print(); }
181  };
182  StatCount statCount;
183 
184 #else
185  struct StatCount {
186  void begin(int) {}
187  void pre(int) {}
188  void end(int) {}
189  void start() {}
190  void noOverlap() {}
191  void overlap() {}
192  void de(float) {}
193  void dp(float) {}
194  };
195  CMS_THREAD_SAFE StatCount statCount;
196 #endif
197 
198 } // namespace
199 
200 namespace {
201  edm::ProductID clusterProductB(const TrackingRecHit* hit) {
202  return reinterpret_cast<const BaseTrackerRecHit*>(hit)->firstClusterRef().id();
203  }
204 } // namespace
205 
207  copyExtras_ = conf.getUntrackedParameter<bool>("copyExtras", true);
208  priorityName_ = conf.getParameter<std::string>("trackAlgoPriorityOrder");
209 
210  std::vector<edm::InputTag> trackProducerTags(conf.getParameter<std::vector<edm::InputTag>>("TrackProducers"));
211  //which of these do I need to turn into vectors?
212  maxNormalizedChisq_ = conf.getParameter<double>("MaxNormalizedChisq");
213  minPT_ = conf.getParameter<double>("MinPT");
214  minFound_ = (unsigned int)conf.getParameter<int>("MinFound");
215  epsilon_ = conf.getParameter<double>("Epsilon");
216  shareFrac_ = conf.getParameter<double>("ShareFrac");
217  allowFirstHitShare_ = conf.getParameter<bool>("allowFirstHitShare");
218  foundHitBonus_ = conf.getParameter<double>("FoundHitBonus");
219  lostHitPenalty_ = conf.getParameter<double>("LostHitPenalty");
220  indivShareFrac_ = conf.getParameter<std::vector<double>>("indivShareFrac");
221  std::string qualityStr = conf.getParameter<std::string>("newQuality");
222  priorityToken = esConsumes<TrackAlgoPriorityOrder, CkfComponentsRecord>(edm::ESInputTag("", priorityName_));
223 
224  if (!qualityStr.empty()) {
226  } else
228 
229  use_sharesInput_ = true;
230  if (epsilon_ > 0.0)
231  use_sharesInput_ = false;
232 
234 
235  for (unsigned int i = 0; i < setsToMerge.size(); i++) {
236  listsToMerge_.push_back(setsToMerge[i].getParameter<std::vector<int>>("tLists"));
237  promoteQuality_.push_back(setsToMerge[i].getParameter<bool>("pQual"));
238  }
239  hasSelector_ = conf.getParameter<std::vector<int>>("hasSelector");
240  copyMVA_ = conf.getParameter<bool>("copyMVA");
241 
242  std::vector<edm::InputTag> selectors(conf.getParameter<std::vector<edm::InputTag>>("selectedTrackQuals"));
243  std::vector<edm::InputTag> mvaStores;
244  if (conf.exists("mvaValueTags")) {
245  mvaStores = conf.getParameter<std::vector<edm::InputTag>>("mvaValueTags");
246  } else {
247  for (int i = 0; i < (int)selectors.size(); i++) {
248  edm::InputTag ntag(selectors[i].label(), "MVAVals");
249  mvaStores.push_back(ntag);
250  }
251  }
252  unsigned int numTrkColl = trackProducerTags.size();
253  if (numTrkColl != hasSelector_.size() || numTrkColl != selectors.size()) {
254  throw cms::Exception("Inconsistent size") << "need same number of track collections and selectors";
255  }
256  if (numTrkColl != hasSelector_.size() || numTrkColl != mvaStores.size()) {
257  throw cms::Exception("Inconsistent size") << "need same number of track collections and MVA stores";
258  }
259  for (unsigned int i = indivShareFrac_.size(); i < numTrkColl; i++) {
260  // edm::LogWarning("TrackListMerger") << "No indivShareFrac for " << trackProducersTags <<". Using default value of 1";
261  indivShareFrac_.push_back(1.0);
262  }
263 
264  trkQualMod_ = conf.getParameter<bool>("writeOnlyTrkQuals");
265  if (trkQualMod_) {
266  bool ok = true;
267  for (unsigned int i = 1; i < numTrkColl; i++) {
268  if (!(trackProducerTags[i] == trackProducerTags[0]))
269  ok = false;
270  }
271  if (!ok) {
272  throw cms::Exception("Bad input") << "to use writeOnlyTrkQuals=True all input InputTags must be the same";
273  }
274  produces<edm::ValueMap<int>>();
275  produces<QualityMaskCollection>("QualityMasks");
276  } else {
277  produces<reco::TrackCollection>();
278 
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  produces<std::vector<Trajectory>>();
290  produces<TrajTrackAssociationCollection>();
291  }
292  produces<edm::ValueMap<float>>("MVAVals");
293  produces<MVACollection>("MVAValues");
294 
295  // Do all the consumes
296  trackProducers_.resize(numTrkColl);
297  for (unsigned int i = 0; i < numTrkColl; ++i) {
298  trackProducers_[i] = hasSelector_[i] > 0 ? edTokens(trackProducerTags[i], selectors[i], mvaStores[i])
299  : edTokens(trackProducerTags[i], mvaStores[i]);
300  }
301 }
302 
303 // Virtual destructor needed.
305 
306 // Functions that gets called by framework every event
308  // extract tracker geometry
309  //
310  //edm::ESHandle<TrackerGeometry> theG;
311  //es.get<TrackerDigiGeometryRecord>().get(theG);
312 
313  // using namespace reco;
314 
316  auto const& trackAlgoPriorityOrder = *priorityH;
317 
318  // get Inputs
319  // if 1 input list doesn't exist, make an empty list, issue a warning, and continue
320  // this allows TrackListMerger to be used as a cleaner only if handed just one list
321  // if both input lists don't exist, will issue 2 warnings and generate an empty output collection
322  //
323  static const reco::TrackCollection s_empty;
324 
325  std::vector<const reco::TrackCollection*> trackColls;
326  std::vector<edm::Handle<reco::TrackCollection>> trackHandles(trackProducers_.size());
327  for (unsigned int i = 0; i < trackProducers_.size(); i++) {
328  trackColls.push_back(nullptr);
329  //edm::Handle<reco::TrackCollection> trackColl;
330  e.getByToken(trackProducers_[i].tk, trackHandles[i]);
331  if (trackHandles[i].isValid()) {
332  trackColls[i] = trackHandles[i].product();
333  } else {
334  edm::LogWarning("TrackListMerger") << "TrackCollection " << trackProducers_[i].tag << " not found";
335  trackColls[i] = &s_empty;
336  }
337  }
338 
339  unsigned int collsSize = trackColls.size();
340  unsigned int rSize = 0;
341  unsigned int trackCollSizes[collsSize];
342  unsigned int trackCollFirsts[collsSize];
343  for (unsigned int i = 0; i != collsSize; i++) {
344  trackCollSizes[i] = trackColls[i]->size();
345  trackCollFirsts[i] = rSize;
346  rSize += trackCollSizes[i];
347  }
348 
349  statCount.begin(rSize);
350 
351  //
352  // quality cuts first
353  //
354  int i = -1;
355 
356  int selected[rSize];
357  int indexG[rSize];
358  bool trkUpdated[rSize];
359  int trackCollNum[rSize];
360  int trackQuals[rSize];
361  float trackMVAs[rSize];
362  reco::TrackBase::TrackAlgorithm oriAlgo[rSize];
363  std::vector<reco::TrackBase::AlgoMask> algoMask(rSize);
364  for (unsigned int j = 0; j < rSize; j++) {
365  indexG[j] = -1;
366  selected[j] = 1;
367  trkUpdated[j] = false;
368  trackCollNum[j] = 0;
369  trackQuals[j] = 0;
370  trackMVAs[j] = -998.0;
372  }
373 
374  int ngood = 0;
375  for (unsigned int j = 0; j != collsSize; j++) {
376  const reco::TrackCollection* tC1 = trackColls[j];
377 
378  edm::Handle<edm::ValueMap<int>> trackSelColl;
379  edm::Handle<edm::ValueMap<float>> trackMVAStore;
380  if (copyMVA_)
381  e.getByToken(trackProducers_[j].tmva, trackMVAStore);
382  if (hasSelector_[j] > 0) {
383  e.getByToken(trackProducers_[j].tsel, trackSelColl);
384  }
385 
386  if (!tC1->empty()) {
387  unsigned int iC = 0;
388  for (reco::TrackCollection::const_iterator track = tC1->begin(); track != tC1->end(); track++) {
389  i++;
390  trackCollNum[i] = j;
391  trackQuals[i] = track->qualityMask();
392  oriAlgo[i] = track->originalAlgo();
393  algoMask[i] = track->algoMask();
394 
395  reco::TrackRef trkRef = reco::TrackRef(trackHandles[j], iC);
396  if (copyMVA_)
397  if ((*trackMVAStore).contains(trkRef.id()))
398  trackMVAs[i] = (*trackMVAStore)[trkRef];
399  if (hasSelector_[j] > 0) {
400  int qual = (*trackSelColl)[trkRef];
401  if (qual < 0) {
402  selected[i] = 0;
403  iC++;
404  continue;
405  } else {
406  trackQuals[i] = qual;
407  }
408  }
409  iC++;
410  selected[i] = trackQuals[i] + 10; //10 is magic number used throughout...
411  if ((short unsigned)track->ndof() < 1) {
412  selected[i] = 0;
413  continue;
414  }
415  if (track->normalizedChi2() > maxNormalizedChisq_) {
416  selected[i] = 0;
417  continue;
418  }
419  if (track->found() < minFound_) {
420  selected[i] = 0;
421  continue;
422  }
423  if (track->pt() < minPT_) {
424  selected[i] = 0;
425  continue;
426  }
427  // good!
428  indexG[i] = ngood++;
429  //if ( beVerb) std::cout << "inverb " << track->pt() << " " << selected[i] << std::endl;
430  } //end loop over tracks
431  } //end more than 0 track
432  } // loop over trackcolls
433 
434  statCount.pre(ngood);
435 
436  //cache the id and rechits of valid hits
437  typedef std::pair<unsigned int, const TrackingRecHit*> IHit;
438  std::vector<std::vector<IHit>> rh1(ngood); // "not an array" of vectors!
439  //const TrackingRecHit* fh1[ngood]; // first hit...
441  float score[ngood];
442 
443  for (unsigned int j = 0; j < rSize; j++) {
444  if (selected[j] == 0)
445  continue;
446  int i = indexG[j];
447  assert(i >= 0);
448  unsigned int collNum = trackCollNum[j];
449  unsigned int trackNum = j - trackCollFirsts[collNum];
450  const reco::Track* track = &((trackColls[collNum])->at(trackNum));
451 
452  algo[i] = track->algo();
453  int validHits = track->numberOfValidHits();
454  int validPixelHits = track->hitPattern().numberOfValidPixelHits();
455  int lostHits = track->numberOfLostHits();
456  score[i] =
457  foundHitBonus_ * validPixelHits + foundHitBonus_ * validHits - lostHitPenalty_ * lostHits - track->chi2();
458 
459  rh1[i].reserve(validHits);
460  auto compById = [](IHit const& h1, IHit const& h2) { return h1.first < h2.first; };
461  for (trackingRecHit_iterator it = track->recHitsBegin(); it != track->recHitsEnd(); ++it) {
462  const TrackingRecHit* hit = (*it);
463  unsigned int id = hit->rawId();
464  if (hit->geographicalId().subdetId() > 2)
465  id &= (~3); // mask mono/stereo in strips...
466  if LIKELY (hit->isValid()) {
467  rh1[i].emplace_back(id, hit);
468  std::push_heap(rh1[i].begin(), rh1[i].end(), compById);
469  }
470  }
471  std::sort_heap(rh1[i].begin(), rh1[i].end(), compById);
472  }
473 
474  //DL here
475  if 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 UNLIKELY (allowFirstHitShare_ && rh1[k1][0].first == rh1[k2][0].first) {
539  const TrackingRecHit* it = rh1[k1][0].second;
540  const TrackingRecHit* jt = rh1[k2][0].second;
541  if (share(it, jt, epsilon_))
542  firstoverlap = 1;
543  }
544 
545  // exploit sorting
546  unsigned int jh = 0;
547  unsigned int ih = 0;
548  while (ih != nh1 && jh != nh2) {
549  // break if not enough to go...
550  // if ( nprecut-noverlap+firstoverlap > int(nh1-ih)) break;
551  // if ( nprecut-noverlap+firstoverlap > int(nh2-jh)) break;
552  auto const id1 = rh1[k1][ih].first;
553  auto const id2 = rh1[k2][jh].first;
554  if (id1 < id2)
555  ++ih;
556  else if (id2 < id1)
557  ++jh;
558  else {
559  // in case of split-hit do full conbinatorics
560  auto li = ih;
561  while ((++li) != nh1 && id1 == rh1[k1][li].first) {
562  }
563  auto lj = jh;
564  while ((++lj) != nh2 && id2 == rh1[k2][lj].first) {
565  }
566  for (auto ii = ih; ii != li; ++ii)
567  for (auto jj = jh; jj != lj; ++jj) {
568  const TrackingRecHit* it = rh1[k1][ii].second;
569  const TrackingRecHit* jt = rh1[k2][jj].second;
570  if (share(it, jt, epsilon_))
571  noverlap++;
572  }
573  jh = lj;
574  ih = li;
575  } // equal ids
576 
577  } //loop over ih & jh
578 
579  bool dupfound =
580  (collNum != collNum2)
581  ? (noverlap - firstoverlap) > (std::min(nhit1, nhit2) - firstoverlap) * shareFrac_
582  : (noverlap - firstoverlap) > (std::min(nhit1, nhit2) - firstoverlap) * indivShareFrac_[collNum];
583 
584  auto seti = [&](unsigned int ii, unsigned int jj) {
585  selected[jj] = 0;
586  selected[ii] = 10 + newQualityMask; // add 10 to avoid the case where mask = 1
587  trkUpdated[ii] = true;
588  if (trackAlgoPriorityOrder.priority(oriAlgo[jj]) < trackAlgoPriorityOrder.priority(oriAlgo[ii]))
589  oriAlgo[ii] = oriAlgo[jj];
590  algoMask[ii] |= algoMask[jj];
591  algoMask[jj] = algoMask[ii]; // in case we keep discarded
592  };
593 
594  if (dupfound) {
595  float score2 = score[k2];
596  constexpr float almostSame =
597  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>>();
648  edm::ValueMap<int>::Filler filler(*vm);
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);
706 
707  if (copyExtras_) {
708  outputTrkExtras = std::make_unique<reco::TrackExtraCollection>();
709  outputTrkExtras->reserve(nToWrite);
711  outputTrkHits = std::make_unique<TrackingRecHitCollection>();
712  outputTrkHits->reserve(nToWrite * 25);
714  if (makeReKeyedSeeds_) {
715  outputSeeds = std::make_unique<TrajectorySeedCollection>();
716  outputSeeds->reserve(nToWrite);
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  TrackingRecHit const& hit = *origSeedRef->recHits().begin();
759  if (hit.isValid()) {
760  edm::ProductID pID = clusterProductB(&hit);
761  // the cluster collection either produced a removalInfo or mot
762  //get the clusterremoval info from the provenance: will rekey if this is found
764  edm::StableProvenance const& prov = e.getStableProvenance(pID);
765  clusterRemovalInfos = edm::InputTag(prov.moduleLabel(), prov.productInstanceName(), prov.processName());
766  doRekeyOnThisSeed = e.getByLabel(clusterRemovalInfos, CRIh);
767  } //valid hit
768  } //nhit!=0
769 
770  if (doRekeyOnThisSeed && !(clusterRemovalInfos == edm::InputTag(""))) {
771  ClusterRemovalRefSetter refSetter(e, clusterRemovalInfos);
772  TrajectorySeed::RecHitContainer newRecHitContainer;
773  newRecHitContainer.reserve(origSeedRef->nHits());
774  for (auto const& recHit : origSeedRef->recHits()) {
775  newRecHitContainer.push_back(recHit);
776  refSetter.reKey(&newRecHitContainer.back());
777  }
778  outputSeeds->push_back(
779  TrajectorySeed(origSeedRef->startingState(), newRecHitContainer, origSeedRef->direction()));
780  }
781  //doRekeyOnThisSeed=true
782  else {
783  //just copy the one we had before
784  outputSeeds->push_back(TrajectorySeed(*origSeedRef));
785  }
787  origSeedRef = edm::RefToBase<TrajectorySeed>(pureRef);
788  } //creating a new seed and rekeying it rechit clusters.
789 
790  // Fill TrackExtra collection
791  outputTrkExtras->push_back(reco::TrackExtra(track->outerPosition(),
792  track->outerMomentum(),
793  track->outerOk(),
794  track->innerPosition(),
795  track->innerMomentum(),
796  track->innerOk(),
797  track->outerStateCovariance(),
798  track->outerDetId(),
799  track->innerStateCovariance(),
800  track->innerDetId(),
801  track->seedDirection(),
802  origSeedRef));
803  seedsRefs[i] = origSeedRef;
804  outputTrks->back().setExtra(reco::TrackExtraRef(refTrkExtras, outputTrkExtras->size() - 1));
805  reco::TrackExtra& tx = outputTrkExtras->back();
806  tx.setResiduals(track->residuals());
807 
808  // fill TrackingRecHits
809  unsigned nh1 = track->recHitsSize();
810  tx.setHits(refTrkHits, outputTrkHits->size(), nh1);
811  tx.setTrajParams(track->extra()->trajParams(), track->extra()->chi2sX5());
812  assert(tx.trajParams().size() == tx.recHitsSize());
813  for (auto hh = track->recHitsBegin(), eh = track->recHitsEnd(); hh != eh; ++hh) {
814  outputTrkHits->push_back((*hh)->clone());
815  }
816  }
817  trackRefs[i] = reco::TrackRef(refTrks, outputTrks->size() - 1);
818 
819  } //end faux loop over tracks
820 
821  //Fill the trajectories, etc. for 1st collection
822  refTrajs = e.getRefBeforePut<std::vector<Trajectory>>();
823 
824  outputTTAss = std::make_unique<TrajTrackAssociationCollection>(refTrajs, refTrks);
825 
826  for (unsigned int ti = 0; ti < trackColls.size(); ti++) {
829  e.getByToken(trackProducers_[ti].traj, hTraj1);
830  e.getByToken(trackProducers_[ti].tass, hTTAss1);
831 
832  if (hTraj1.failedToGet() || hTTAss1.failedToGet())
833  continue;
834 
835  for (size_t i = 0, n = hTraj1->size(); i < n; ++i) {
836  edm::Ref<std::vector<Trajectory>> trajRef(hTraj1, i);
837  TrajTrackAssociationCollection::const_iterator match = hTTAss1->find(trajRef);
838  if (match != hTTAss1->end()) {
839  const edm::Ref<reco::TrackCollection>& trkRef = match->val;
840  uint32_t oldKey = trackCollFirsts[ti] + static_cast<uint32_t>(trkRef.key());
841  if (trackRefs[oldKey].isNonnull()) {
842  outputTrajs->push_back(*trajRef);
843  //if making extras and the seeds at the same time, change the seed ref on the trajectory
845  outputTrajs->back().setSeedRef(seedsRefs[oldKey]);
846  outputTTAss->insert(edm::Ref<std::vector<Trajectory>>(refTrajs, outputTrajs->size() - 1), trackRefs[oldKey]);
847  }
848  }
849  }
850  }
851 
852  statCount.end(outputTrks->size());
853 
854  edm::ProductID nPID = refTrks.id();
855  edm::TestHandle<reco::TrackCollection> outHandle(outputTrks.get(), nPID);
856  fillerMVA.insert(outHandle, mvaVec.begin(), mvaVec.end());
857  fillerMVA.fill();
858 
860  if (copyMVA_) {
861  e.put(std::move(vmMVA), "MVAVals");
862  auto mvas = std::make_unique<MVACollection>(mvaVec.begin(), mvaVec.end());
863  e.put(std::move(mvas), "MVAValues");
864  }
865  if (copyExtras_) {
868  if (makeReKeyedSeeds_)
870  }
872  e.put(std::move(outputTTAss));
873  return;
874 
875 } //end produce
876 
879 
PropagationDirection direction() const
const edm::RefToBase< TrajectorySeed > & seedRef() const
Definition: Track.h:155
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:133
std::vector< bool > promoteQuality_
const edm::EventSetup & c
std::unique_ptr< TrajTrackAssociationCollection > outputTTAss
const bool isValid(const Frame &aFrame, const FrameQuality &aQuality, const uint16_t aExpectedPos)
WorkSpace int float eps
const TrackExtraRef & extra() const
reference to &quot;extra&quot; object
Definition: Track.h:139
reco::TrackRefProd refTrks
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
TkEDGetTokenss edTokens(const edm::InputTag &tag, const edm::InputTag &mvatag)
TrackQuality
track quality
Definition: TrackBase.h:150
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
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:34
void reKey(TrackingRecHit *hit) const
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:801
key_type key() const
Accessor for product key.
Definition: Ref.h:250
~TrackListMerger() override
int ii
Definition: cuy.py:589
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
assert(be >=bs)
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_
bool overlap(const reco::Muon &muon1, const reco::Muon &muon2, double pullX=1.0, double pullY=1.0, bool checkAdjacentChambers=false)
const math::XYZPoint & innerPosition() const
position of the innermost hit
Definition: Track.h:56
TrackAlgorithm algo() const
Definition: TrackBase.h:547
tuple d
Definition: ztail.py:151
char const * label
void push_back(D *&d)
Definition: OwnVector.h:326
std::string const & productInstanceName() const
std::unique_ptr< TrackingRecHitCollection > outputTrkHits
edm::EDGetTokenT< edm::ValueMap< int > > tsel
std::vector< TrajectorySeed > TrajectorySeedCollection
double chi2() const
chi-squared of the fit
Definition: TrackBase.h:587
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
def move
Definition: eostools.py:511
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
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
TkEDGetTokenss edTokens(const edm::InputTag &tag, const edm::InputTag &seltag, const edm::InputTag &mvatag)
#define CMS_THREAD_SAFE
edm::ESGetToken< TrackAlgoPriorityOrder, CkfComponentsRecord > priorityToken
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:798
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
edm::RefProd< std::vector< Trajectory > > refTrajs
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:500
RefProd< PROD > getRefBeforePut()
Definition: Event.h:158
std::unique_ptr< reco::TrackExtraCollection > outputTrkExtras
std::vector< int > hasSelector_
#define dso_hidden
Definition: Visibility.h:12
bool failedToGet() const
Definition: HandleBase.h:72
RecHitRange recHits() const
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:504
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
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
bool isValid() const
constexpr std::array< unsigned int, reco::TrackBase::algoSize > trackAlgoPriorityOrder
unsigned int nHits() const
unsigned int minFound_
std::vector< TkEDGetTokenss > trackProducers_
const math::XYZVector & innerMomentum() const
momentum vector at the innermost hit position
Definition: Track.h:59
std::string const & processName() const
ProductID id() const
Accessor for product ID.
Definition: RefProd.h:124
string end
Definition: dataset.py:937
const TrackResiduals & residuals() const
get the residuals
Definition: Track.h:158
T begin() const
Definition: Range.h:15
int numberOfValidPixelHits() const
Definition: HitPattern.h:825
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
tuple cout
Definition: gather_cfg.py:144
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
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:157
Log< level::Warning, false > LogWarning
std::string const & moduleLabel() const
std::string priorityName_
StableProvenance const & getStableProvenance(BranchID const &theID) const
Definition: Event.cc:124
T x() const
Definition: PV3DBase.h:59
virtual LocalPoint localPosition() const =0
unsigned int innerDetId() const
DetId of the detector on which surface the innermost state is located.
Definition: Track.h:82
id_type rawId() const
auto const & hh
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