CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrackListMerger.cc
Go to the documentation of this file.
1 //
2 // Package: RecoTracker/FinalTrackSelectors
3 // Class: TrackListMerger
4 //
5 // Description: TrackList Cleaner and Merger able to deal with many lists
6 //
7 // Original Author: David Lange
8 // Created: April 4, 2011
9 //
10 
11 
12 #include <memory>
13 #include <string>
14 #include <iostream>
15 #include <cmath>
16 #include <vector>
17 
19 
27 
29 
34 
35 //#include "DataFormats/TrackReco/src/classes.h"
36 
38 
39 #ifdef STAT_TSB
40 #include <x86intrin.h>
41 #endif
42 
43 namespace {
44 #ifdef STAT_TSB
45 inline volatile unsigned long long rdtsc() {
46  return __rdtsc();
47 }
48 
49  struct StatCount {
50  float maxDP=0.;
51  float maxDE=0.;
52  unsigned long long st;
53  long long totBegin=0;
54  long long totPre=0;
55  long long totEnd=0;
56  unsigned long long timeNo; // no-overlap
57  unsigned long long timeOv; // overlap
58  void begin(int tt) {
59  totBegin+=tt;
60  }
61  void start() { st=rdtsc(); }
62  void noOverlap() { timeNo += (rdtsc()-st);}
63  void overlap() { timeOv += (rdtsc()-st);}
64  void pre(int tt) { totPre+=tt;}
65  void end(int tt) { totEnd+=tt;}
66  void de(float d) { if (d>maxDE) maxDE=d;}
67  void dp(float d) { if (d>maxDP) maxDP=d;}
68 
69 
70  void print() const {
71  std::cout << "TrackListMerger stat\nBegin/Pre/End/maxDPhi/maxDEta/Overlap/NoOverlap "
72  << totBegin <<'/'<< totPre <<'/'<< totEnd <<'/'<< maxDP <<'/'<< maxDE
73  <<'/'<< timeOv/1000 <<'/'<< timeNo/1000
74  << std::endl;
75  }
76  StatCount() {}
77  ~StatCount() { print();}
78  };
79  StatCount statCount;
80 
81 #else
82  struct StatCount {
83  void begin(int){}
84  void pre(int){}
85  void end(int){}
86  void start(){}
87  void noOverlap(){}
88  void overlap(){}
89  void de(float){}
90  void dp(float){}
91 
92 
93  };
94  [[cms::thread_safe]] StatCount statCount;
95 #endif
96 
97 
98 }
99 
100 
101 
102 namespace cms
103 {
104 
106  return reinterpret_cast<const BaseTrackerRecHit *>(hit)->firstClusterRef().id();
107  }
108 
109 
110 
112  copyExtras_ = conf.getUntrackedParameter<bool>("copyExtras", true);
113 
114  std::vector<edm::InputTag> trackProducerTags(conf.getParameter<std::vector<edm::InputTag> >("TrackProducers"));
115  //which of these do I need to turn into vectors?
116  maxNormalizedChisq_ = conf.getParameter<double>("MaxNormalizedChisq");
117  minPT_ = conf.getParameter<double>("MinPT");
118  minFound_ = (unsigned int)conf.getParameter<int>("MinFound");
119  epsilon_ = conf.getParameter<double>("Epsilon");
120  shareFrac_ = conf.getParameter<double>("ShareFrac");
121  allowFirstHitShare_ = conf.getParameter<bool>("allowFirstHitShare");
122  foundHitBonus_ = conf.getParameter<double>("FoundHitBonus");
123  lostHitPenalty_ = conf.getParameter<double>("LostHitPenalty");
124  indivShareFrac_=conf.getParameter<std::vector<double> >("indivShareFrac");
125  std::string qualityStr = conf.getParameter<std::string>("newQuality");
126 
127  if (qualityStr != "") {
129  }
130  else
132 
133  use_sharesInput_ = true;
134  if ( epsilon_ > 0.0 )use_sharesInput_ = false;
135 
136  edm::VParameterSet setsToMerge=conf.getParameter<edm::VParameterSet>("setsToMerge");
137 
138  for ( unsigned int i=0; i<setsToMerge.size(); i++) {
139  listsToMerge_.push_back(setsToMerge[i].getParameter<std::vector< int> >("tLists"));
140  promoteQuality_.push_back(setsToMerge[i].getParameter<bool>("pQual"));
141  }
142  hasSelector_= conf.getParameter<std::vector<int> >("hasSelector");
143  std::vector<edm::InputTag> selectors(conf.getParameter<std::vector<edm::InputTag> >("selectedTrackQuals"));
144  std::vector<edm::InputTag> mvaStores;
145  if(conf.exists("mvaValueTags")){
146  mvaStores = conf.getParameter<std::vector<edm::InputTag> >("mvaValueTags");
147  }else{
148  for (int i = 0; i < (int)selectors.size(); i++){
149  edm::InputTag ntag(selectors[i].label(),"MVAVals");
150  mvaStores.push_back(ntag);
151  }
152  }
153  unsigned int numTrkColl = trackProducerTags.size();
154  if (numTrkColl != hasSelector_.size() || numTrkColl != selectors.size()) {
155  throw cms::Exception("Inconsistent size") << "need same number of track collections and selectors";
156  }
157  if (numTrkColl != hasSelector_.size() || numTrkColl != mvaStores.size()) {
158  throw cms::Exception("Inconsistent size") << "need same number of track collections and MVA stores";
159  }
160  for (unsigned int i = indivShareFrac_.size(); i < numTrkColl; i++) {
161  // edm::LogWarning("TrackListMerger") << "No indivShareFrac for " << trackProducersTags <<". Using default value of 1";
162  indivShareFrac_.push_back(1.0);
163  }
164 
165  trkQualMod_=conf.getParameter<bool>("writeOnlyTrkQuals");
166  if ( trkQualMod_) {
167  bool ok=true;
168  for ( unsigned int i=1; i<numTrkColl; i++) {
169  if (!(trackProducerTags[i]==trackProducerTags[0])) ok=false;
170  }
171  if ( !ok) {
172  throw cms::Exception("Bad input") << "to use writeOnlyTrkQuals=True all input InputTags must be the same";
173  }
174  produces<edm::ValueMap<int> >();
175  }
176  else{
177  produces<reco::TrackCollection>();
178 
179  makeReKeyedSeeds_ = conf.getUntrackedParameter<bool>("makeReKeyedSeeds",false);
180  if (makeReKeyedSeeds_){
181  copyExtras_=true;
182  produces<TrajectorySeedCollection>();
183  }
184 
185  if (copyExtras_) {
186  produces<reco::TrackExtraCollection>();
187  produces<TrackingRecHitCollection>();
188  }
189  produces< std::vector<Trajectory> >();
190  produces< TrajTrackAssociationCollection >();
191  }
192  produces<edm::ValueMap<float> >("MVAVals");
193 
194  // Do all the consumes
195  trackProducers_.resize(numTrkColl);
196  for (unsigned int i = 0; i < numTrkColl; ++i) {
197  trackProducers_[i] = hasSelector_[i]>0 ? edTokens(trackProducerTags[i], selectors[i], mvaStores[i]) : edTokens(trackProducerTags[i], mvaStores[i]);
198  }
199  }
200 
201 
202  // Virtual destructor needed.
204 
205  // Functions that gets called by framework every event
207  {
208  // extract tracker geometry
209  //
210  //edm::ESHandle<TrackerGeometry> theG;
211  //es.get<TrackerDigiGeometryRecord>().get(theG);
212 
213  // using namespace reco;
214 
215  // get Inputs
216  // if 1 input list doesn't exist, make an empty list, issue a warning, and continue
217  // this allows TrackListMerger to be used as a cleaner only if handed just one list
218  // if both input lists don't exist, will issue 2 warnings and generate an empty output collection
219  //
220  static const reco::TrackCollection s_empty;
221 
222  std::vector<const reco::TrackCollection *> trackColls;
223  std::vector<edm::Handle<reco::TrackCollection> > trackHandles(trackProducers_.size());
224  for ( unsigned int i=0; i<trackProducers_.size(); i++) {
225  trackColls.push_back(0);
226  //edm::Handle<reco::TrackCollection> trackColl;
227  e.getByToken(trackProducers_[i].tk, trackHandles[i]);
228  if (trackHandles[i].isValid()) {
229  trackColls[i]= trackHandles[i].product();
230  } else {
231  edm::LogWarning("TrackListMerger") << "TrackCollection " << trackProducers_[i].tag <<" not found";
232  trackColls[i]=&s_empty;
233  }
234  }
235 
236  unsigned int collsSize =trackColls.size();
237  unsigned int rSize=0;
238  unsigned int trackCollSizes[collsSize];
239  unsigned int trackCollFirsts[collsSize];
240  for (unsigned int i=0; i!=collsSize; i++) {
241  trackCollSizes[i]=trackColls[i]->size();
242  trackCollFirsts[i]=rSize;
243  rSize+=trackCollSizes[i];
244  }
245 
246  statCount.begin(rSize);
247 
248  //
249  // quality cuts first
250  //
251  int i=-1;
252 
253  int selected[rSize];
254  int indexG[rSize];
255  bool trkUpdated[rSize];
256  int trackCollNum[rSize];
257  int trackQuals[rSize];
258  float trackMVAs[rSize];
259  for (unsigned int j=0; j<rSize;j++) {
260  indexG[j]=-1; selected[j]=1; trkUpdated[j]=false; trackCollNum[j]=0; trackQuals[j]=0;trackMVAs[j] = -99.0;
261  }
262 
263  int ngood=0;
264  for (unsigned int j=0; j!= collsSize; j++) {
265  const reco::TrackCollection *tC1=trackColls[j];
266 
267  edm::Handle<edm::ValueMap<int> > trackSelColl;
268  edm::Handle<edm::ValueMap<float> > trackMVAStore;
269  e.getByToken(trackProducers_[j].tmva, trackMVAStore);
270  if ( hasSelector_[j]>0 ){
271  e.getByToken(trackProducers_[j].tsel, trackSelColl);
272  }
273 
274  if ( 0<tC1->size() ){
275  unsigned int iC=0;
276  for (reco::TrackCollection::const_iterator track=tC1->begin(); track!=tC1->end(); track++){
277  i++;
278  trackCollNum[i]=j;
279  trackQuals[i]=track->qualityMask();
280 
281  reco::TrackRef trkRef=reco::TrackRef(trackHandles[j],iC);
282  if((*trackMVAStore).contains(trkRef.id()))trackMVAs[i] = (*trackMVAStore)[trkRef];
283  if ( hasSelector_[j]>0 ) {
284  int qual=(*trackSelColl)[trkRef];
285  if ( qual < 0 ) {
286  selected[i]=0;
287  iC++;
288  continue;
289  }else{
290  trackQuals[i]=qual;
291  }
292  }
293  iC++;
294  selected[i]=trackQuals[i]+10;//10 is magic number used throughout...
295  if ((short unsigned)track->ndof() < 1){
296  selected[i]=0;
297  continue;
298  }
299  if (track->normalizedChi2() > maxNormalizedChisq_){
300  selected[i]=0;
301  continue;
302  }
303  if (track->found() < minFound_){
304  selected[i]=0;
305  continue;
306  }
307  if (track->pt() < minPT_){
308  selected[i]=0;
309  continue;
310  }
311  // good!
312  indexG[i] = ngood++;
313  //if ( beVerb) std::cout << "inverb " << track->pt() << " " << selected[i] << std::endl;
314  }//end loop over tracks
315  }//end more than 0 track
316  } // loop over trackcolls
317 
318 
319  statCount.pre(ngood);
320 
321  //cache the id and rechits of valid hits
322  typedef std::pair<unsigned int, const TrackingRecHit*> IHit;
323  std::vector<IHit> rh1[ngood]; // an array of vectors!
324  //const TrackingRecHit* fh1[ngood]; // first hit...
325  uint8_t algo[ngood];
326  float score[ngood];
327 
328 
329  for ( unsigned int j=0; j<rSize; j++) {
330  if (selected[j]==0) continue;
331  int i = indexG[j];
332  assert(i>=0);
333  unsigned int collNum=trackCollNum[j];
334  unsigned int trackNum=j-trackCollFirsts[collNum];
335  const reco::Track *track=&((trackColls[collNum])->at(trackNum));
336 
337  algo[i]=track->algo();
338  int validHits=track->numberOfValidHits();
339  int lostHits=track->numberOfLostHits();
340  score[i] = foundHitBonus_*validHits - lostHitPenalty_*lostHits - track->chi2();
341 
342 
343  rh1[i].reserve(validHits) ;
344  auto compById = [](IHit const & h1, IHit const & h2) {return h1.first < h2.first;};
345  // fh1[i] = &(**track->recHitsBegin());
346  for (trackingRecHit_iterator it = track->recHitsBegin(); it != track->recHitsEnd(); ++it) {
347  const TrackingRecHit* hit = &(**it);
348  unsigned int id = hit->rawId() ;
349  if(hit->geographicalId().subdetId()>2) id &= (~3); // mask mono/stereo in strips...
350  if likely(hit->isValid()) { rh1[i].emplace_back(id,hit); std::push_heap(rh1[i].begin(),rh1[i].end(),compById); }
351  }
352  std::sort_heap(rh1[i].begin(),rh1[i].end(),compById);
353  }
354 
355  //DL here
356  if likely(ngood>1 && collsSize>1)
357  for ( unsigned int ltm=0; ltm<listsToMerge_.size(); ltm++) {
358  int saveSelected[rSize];
359  bool notActive[collsSize];
360  for (unsigned int cn=0;cn!=collsSize;++cn)
361  notActive[cn]= find(listsToMerge_[ltm].begin(),listsToMerge_[ltm].end(),cn)==listsToMerge_[ltm].end();
362 
363  for ( unsigned int i=0; i<rSize; i++) saveSelected[i]=selected[i];
364 
365  //DL protect against 0 tracks?
366  for ( unsigned int i=0; i<rSize-1; i++) {
367  if (selected[i]==0) continue;
368  unsigned int collNum=trackCollNum[i];
369 
370  //check that this track is in one of the lists for this iteration
371  if (notActive[collNum]) continue;
372 
373  int k1 = indexG[i];
374  unsigned int nh1=rh1[k1].size();
375  int qualityMaskT1 = trackQuals[i];
376 
377  int nhit1 = nh1; // validHits[k1];
378  float score1 = score[k1];
379 
380  // start at next collection
381  for ( unsigned int j=i+1; j<rSize; j++) {
382  if (selected[j]==0) continue;
383  unsigned int collNum2=trackCollNum[j];
384  if ( (collNum == collNum2) && indivShareFrac_[collNum] > 0.99) continue;
385  //check that this track is in one of the lists for this iteration
386  if (notActive[collNum2]) continue;
387 
388  int k2 = indexG[j];
389 
390  int newQualityMask = -9; //avoid resetting quality mask if not desired 10+ -9 =1
391  if (promoteQuality_[ltm]) {
392  int maskT1= saveSelected[i]>1? saveSelected[i]-10 : qualityMaskT1;
393  int maskT2= saveSelected[j]>1? saveSelected[j]-10 : trackQuals[j];
394  newQualityMask =(maskT1 | maskT2); // take OR of trackQuality
395  }
396  unsigned int nh2=rh1[k2].size();
397  int nhit2 = nh2;
398 
399 
400  auto share = use_sharesInput_ ?
401  [](const TrackingRecHit* it,const TrackingRecHit* jt, float)->bool { return it->sharesInput(jt,TrackingRecHit::some); } :
402  [](const TrackingRecHit* it,const TrackingRecHit* jt, float eps)->bool {
403  float delta = std::abs ( it->localPosition().x()-jt->localPosition().x() );
404  return (it->geographicalId()==jt->geographicalId())&&(delta<eps);
405  };
406 
407  statCount.start();
408 
409  //loop over rechits
410  int noverlap=0;
411  int firstoverlap=0;
412  // check first hit (should use REAL first hit?)
413  if unlikely(allowFirstHitShare_ && rh1[k1][0].first==rh1[k2][0].first ) {
414  const TrackingRecHit* it = rh1[k1][0].second;
415  const TrackingRecHit* jt = rh1[k2][0].second;
416  if (share(it,jt,epsilon_)) firstoverlap=1;
417  }
418 
419 
420  // exploit sorting
421  unsigned int jh=0;
422  unsigned int ih=0;
423  while (ih!=nh1 && jh!=nh2) {
424  // break if not enough to go...
425  // if ( nprecut-noverlap+firstoverlap > int(nh1-ih)) break;
426  // if ( nprecut-noverlap+firstoverlap > int(nh2-jh)) break;
427  auto const id1 = rh1[k1][ih].first;
428  auto const id2 = rh1[k2][jh].first;
429  if (id1<id2) ++ih;
430  else if (id2<id1) ++jh;
431  else {
432  // in case of split-hit do full conbinatorics
433  auto li=ih; while( (++li)!=nh1 && id1 == rh1[k1][li].first);
434  auto lj=jh; while( (++lj)!=nh2 && id2 == rh1[k2][lj].first);
435  for (auto ii=ih; ii!=li; ++ii)
436  for (auto jj=jh; jj!=lj; ++jj) {
437  const TrackingRecHit* it = rh1[k1][ii].second;
438  const TrackingRecHit* jt = rh1[k2][jj].second;
439  if (share(it,jt,epsilon_)) noverlap++;
440  }
441  jh=lj; ih=li;
442  } // equal ids
443 
444  } //loop over ih & jh
445 
446  bool dupfound = (collNum != collNum2) ? (noverlap-firstoverlap) > (std::min(nhit1,nhit2)-firstoverlap)*shareFrac_ :
447  (noverlap-firstoverlap) > (std::min(nhit1,nhit2)-firstoverlap)*indivShareFrac_[collNum];
448 
449  if ( dupfound ) {
450  float score2 = score[k2];
451  constexpr float almostSame = 0.01f; // difference rather than ratio due to possible negative values for score
452  if ( score1 - score2 > almostSame ) {
453  selected[j]=0;
454  selected[i]=10+newQualityMask; // add 10 to avoid the case where mask = 1
455  trkUpdated[i]=true;
456  } else if ( score2 - score1 > almostSame ) {
457  selected[i]=0;
458  selected[j]=10+newQualityMask; // add 10 to avoid the case where mask = 1
459  trkUpdated[j]=true;
460  }else{
461  // If tracks from both iterations are virtually identical, choose the one with the best quality or with lower algo
464  //same quality, pick earlier algo
465  if (algo[k1] <= algo[k2]) {
466  selected[j]=0;
467  selected[i]=10+newQualityMask; // add 10 to avoid the case where mask = 1
468  trkUpdated[i]=true;
469  } else {
470  selected[i]=0;
471  selected[j]=10+newQualityMask; // add 10 to avoid the case where mask = 1
472  trkUpdated[j]=true;
473  }
474  } else if ((trackQuals[j] & (1<<reco::TrackBase::loose|1<<reco::TrackBase::tight|1<<reco::TrackBase::highPurity) ) <
476  selected[j]=0;
477  selected[i]=10+newQualityMask; // add 10 to avoid the case where mask = 1
478  trkUpdated[i]=true;
479  }else{
480  selected[i]=0;
481  selected[j]=10+newQualityMask; // add 10 to avoid the case where mask = 1
482  trkUpdated[j]=true;
483  }
484  }//end fi < fj
485  statCount.overlap();
486  /*
487  if (at0[k1]&&at0[k2]) {
488  statCount.dp(dphi);
489  if (dz<1.f) statCount.de(deta);
490  }
491  */
492  }//end got a duplicate
493  else {
494  statCount.noOverlap();
495  }
496  //stop if the ith track is now unselected
497  if (selected[i]==0) break;
498  }//end track2 loop
499  }//end track loop
500  } //end loop over track list sets
501 
502 
503 
504  std::auto_ptr<edm::ValueMap<float> > vmMVA = std::auto_ptr<edm::ValueMap<float> >(new edm::ValueMap<float>);
505  edm::ValueMap<float>::Filler fillerMVA(*vmMVA);
506 
507 
508 
509  // special case - if just doing the trkquals
510  if (trkQualMod_) {
511  unsigned int tSize=trackColls[0]->size();
512  std::auto_ptr<edm::ValueMap<int> > vm = std::auto_ptr<edm::ValueMap<int> >(new edm::ValueMap<int>);
513  edm::ValueMap<int>::Filler filler(*vm);
514 
515  std::vector<int> finalQuals(tSize,-1); //default is unselected
516  for ( unsigned int i=0; i<rSize; i++) {
517  unsigned int tNum=i%tSize;
518 
519  if (selected[i]>1 ) {
520  finalQuals[tNum]=selected[i]-10;
521  if (trkUpdated[i])
522  finalQuals[tNum]=(finalQuals[tNum] | (1<<qualityToSet_));
523  }
524  if ( selected[i]==1 )
525  finalQuals[tNum]=trackQuals[i];
526  }
527 
528  filler.insert(trackHandles[0], finalQuals.begin(),finalQuals.end());
529  filler.fill();
530 
531  e.put(vm);
532 
533  std::vector<float> mvaVec(tSize,-99);
534 
535  for(unsigned int i = 0; i < rSize; i++){
536  unsigned int tNum = i % tSize;
537  mvaVec[tNum] = trackMVAs[tNum];
538  }
539 
540  fillerMVA.insert(trackHandles[0],mvaVec.begin(),mvaVec.end());
541  fillerMVA.fill();
542  e.put(vmMVA,"MVAVals");
543  return;
544  }
545 
546 
547  //
548  // output selected tracks - if any
549  //
550 
551  reco::TrackRef trackRefs[rSize];
552  edm::RefToBase<TrajectorySeed> seedsRefs[rSize];
553 
554  unsigned int nToWrite=0;
555  for ( unsigned int i=0; i<rSize; i++)
556  if (selected[i]!=0) nToWrite++;
557 
558  std::vector<float> mvaVec;
559 
560  outputTrks = std::auto_ptr<reco::TrackCollection>(new reco::TrackCollection);
561  outputTrks->reserve(nToWrite);
563 
564  if (copyExtras_) {
565  outputTrkExtras = std::auto_ptr<reco::TrackExtraCollection>(new reco::TrackExtraCollection);
566  outputTrkExtras->reserve(nToWrite);
568  outputTrkHits = std::auto_ptr<TrackingRecHitCollection>(new TrackingRecHitCollection);
569  outputTrkHits->reserve(nToWrite*25);
571  if (makeReKeyedSeeds_){
572  outputSeeds = std::auto_ptr<TrajectorySeedCollection>(new TrajectorySeedCollection);
573  outputSeeds->reserve(nToWrite);
575  }
576  }
577 
578 
579  outputTrajs = std::auto_ptr< std::vector<Trajectory> >(new std::vector<Trajectory>());
580  outputTrajs->reserve(rSize);
581  outputTTAss = std::auto_ptr< TrajTrackAssociationCollection >(new TrajTrackAssociationCollection());
582 
583 
584 
585  for ( unsigned int i=0; i<rSize; i++) {
586  if (selected[i]==0) {
587  trackRefs[i]=reco::TrackRef();
588  continue;
589  }
590 
591  unsigned int collNum=trackCollNum[i];
592  unsigned int trackNum=i-trackCollFirsts[collNum];
593  const reco::Track *track=&((trackColls[collNum])->at(trackNum));
594  outputTrks->push_back( reco::Track( *track ) );
595  mvaVec.push_back(trackMVAs[i]);
596  if (selected[i]>1 ) {
597  outputTrks->back().setQualityMask(selected[i]-10);
598  if (trkUpdated[i])
599  outputTrks->back().setQuality(qualityToSet_);
600  }
601  //might duplicate things, but doesnt hurt
602  if ( selected[i]==1 )
603  outputTrks->back().setQualityMask(trackQuals[i]);
604 
605  // if ( beVerb ) std::cout << "selected " << outputTrks->back().pt() << " " << outputTrks->back().qualityMask() << " " << selected[i] << std::endl;
606 
607  //fill the TrackCollection
608  if (copyExtras_) {
609  edm::RefToBase<TrajectorySeed> origSeedRef = track->seedRef();
610  //creating a seed with rekeyed clusters if required
611  if (makeReKeyedSeeds_){
612  bool doRekeyOnThisSeed=false;
613 
614  edm::InputTag clusterRemovalInfos("");
615  //grab on of the hits of the seed
616  if (origSeedRef->nHits()!=0){
617  TrajectorySeed::const_iterator firstHit=origSeedRef->recHits().first;
618  const TrackingRecHit *hit = &*firstHit;
619  if (firstHit->isValid()){
621  // the cluster collection either produced a removalInfo or mot
622  //get the clusterremoval info from the provenance: will rekey if this is found
624  edm::Provenance prov=e.getProvenance(pID);
625  clusterRemovalInfos=edm::InputTag(prov.moduleLabel(),
626  prov.productInstanceName(),
627  prov.processName());
628  doRekeyOnThisSeed=e.getByLabel(clusterRemovalInfos,CRIh);
629  }//valid hit
630  }//nhit!=0
631 
632  if (doRekeyOnThisSeed && !(clusterRemovalInfos==edm::InputTag(""))) {
633  ClusterRemovalRefSetter refSetter(e,clusterRemovalInfos);
634  TrajectorySeed::recHitContainer newRecHitContainer;
635  newRecHitContainer.reserve(origSeedRef->nHits());
636  TrajectorySeed::const_iterator iH=origSeedRef->recHits().first;
637  TrajectorySeed::const_iterator iH_end=origSeedRef->recHits().second;
638  for (;iH!=iH_end;++iH){
639  newRecHitContainer.push_back(*iH);
640  refSetter.reKey(&newRecHitContainer.back());
641  }
642  outputSeeds->push_back( TrajectorySeed( origSeedRef->startingState(),
643  newRecHitContainer,
644  origSeedRef->direction()));
645  }
646  //doRekeyOnThisSeed=true
647  else{
648  //just copy the one we had before
649  outputSeeds->push_back( TrajectorySeed(*origSeedRef));
650  }
652  origSeedRef=edm::RefToBase<TrajectorySeed>( pureRef);
653  }//creating a new seed and rekeying it rechit clusters.
654 
655  // Fill TrackExtra collection
656  outputTrkExtras->push_back( reco::TrackExtra(
657  track->outerPosition(), track->outerMomentum(), track->outerOk(),
658  track->innerPosition(), track->innerMomentum(), track->innerOk(),
659  track->outerStateCovariance(), track->outerDetId(),
660  track->innerStateCovariance(), track->innerDetId(),
661  track->seedDirection(), origSeedRef ) );
662  seedsRefs[i]=origSeedRef;
663  outputTrks->back().setExtra( reco::TrackExtraRef( refTrkExtras, outputTrkExtras->size() - 1) );
664  reco::TrackExtra & tx = outputTrkExtras->back();
665  tx.setResiduals(track->residuals());
666 
667  // fill TrackingRecHits
668  unsigned nh1=track->recHitsSize();
669  for ( unsigned ih=0; ih<nh1; ++ih ) {
670  //const TrackingRecHit*hit=&((*(track->recHit(ih))));
671  outputTrkHits->push_back( track->recHit(ih)->clone() );
672  tx.add( TrackingRecHitRef( refTrkHits, outputTrkHits->size() - 1) );
673  }
674  }
675  trackRefs[i] = reco::TrackRef(refTrks, outputTrks->size() - 1);
676 
677 
678  }//end faux loop over tracks
679 
680  //Fill the trajectories, etc. for 1st collection
681  refTrajs = e.getRefBeforePut< std::vector<Trajectory> >();
682 
683  for (unsigned int ti=0; ti<trackColls.size(); ti++) {
686  e.getByToken(trackProducers_[ti].traj, hTraj1);
687  e.getByToken(trackProducers_[ti].tass, hTTAss1);
688 
689  if (hTraj1.failedToGet() || hTTAss1.failedToGet()) continue;
690 
691  for (size_t i = 0, n = hTraj1->size(); i < n; ++i) {
692  edm::Ref< std::vector<Trajectory> > trajRef(hTraj1, i);
693  TrajTrackAssociationCollection::const_iterator match = hTTAss1->find(trajRef);
694  if (match != hTTAss1->end()) {
695  const edm::Ref<reco::TrackCollection> &trkRef = match->val;
696  short oldKey = trackCollFirsts[ti]+static_cast<short>(trkRef.key());
697  if (trackRefs[oldKey].isNonnull()) {
698  outputTrajs->push_back( *trajRef );
699  //if making extras and the seeds at the same time, change the seed ref on the trajectory
701  outputTrajs->back().setSeedRef( seedsRefs[oldKey] );
702  outputTTAss->insert ( edm::Ref< std::vector<Trajectory> >(refTrajs, outputTrajs->size() - 1),
703  trackRefs[oldKey] );
704  }
705  }
706  }
707  }
708 
709  statCount.end(outputTrks->size());
710 
711  edm::ProductID nPID = refTrks.id();
713  fillerMVA.insert(outHandle,mvaVec.begin(),mvaVec.end());
714  fillerMVA.fill();
715 
716  e.put(outputTrks);
717  e.put(vmMVA,"MVAVals");
718  if (copyExtras_) {
719  e.put(outputTrkExtras);
720  e.put(outputTrkHits);
721  if (makeReKeyedSeeds_)
722  e.put(outputSeeds);
723  }
724  e.put(outputTrajs);
725  e.put(outputTTAss);
726  return;
727 
728  }//end produce
729 }
PropagationDirection direction() const
dbl * delta
Definition: mlp_gen.cc:36
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
virtual void produce(edm::Event &e, const edm::EventSetup &c) override
reference back()
Definition: OwnVector.h:321
tuple start
Check for commandline option errors.
Definition: dqm_diff.py:58
std::auto_ptr< TrajectorySeedCollection > outputSeeds
std::auto_ptr< std::vector< Trajectory > > outputTrajs
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
size_t recHitsSize() const
Get number of RecHits. (Warning, this includes invalid hits, which are not physical hits)...
Definition: Track.h:69
TkEDGetTokenss edTokens(const edm::InputTag &tag, const edm::InputTag &seltag, const edm::InputTag &mvatag)
std::vector< ParameterSet > VParameterSet
Definition: ParameterSet.h:33
std::string print(const Track &, edm::Verbosity=edm::Concise)
Track print utility.
Definition: print.cc:8
void reKey(TrackingRecHit *hit) const
std::auto_ptr< reco::TrackCollection > outputTrks
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:10
bool exists(std::string const &parameterName) const
checks if a parameter exists
bool innerOk() const
return true if the innermost hit is valid
Definition: Track.h:40
unsigned short numberOfLostHits() const
number of cases where track crossed a layer without getting a hit.
Definition: TrackBase.h:232
std::vector< std::vector< int > > listsToMerge_
int ii
Definition: cuy.py:588
std::string const & processName() const
Definition: Provenance.h:61
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
const math::XYZPoint & outerPosition() const
position of the outermost hit
Definition: Track.h:47
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:42
TrackAlgorithm algo() const
Definition: TrackBase.h:330
#define unlikely(x)
Definition: Likely.h:21
void push_back(D *&d)
Definition: OwnVector.h:274
std::vector< TrajectorySeed > TrajectorySeedCollection
double chi2() const
chi-squared of the fit
Definition: TrackBase.h:105
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
recHitContainer::const_iterator const_iterator
edm::Ref< TrackingRecHitCollection > TrackingRecHitRef
persistent reference to a TrackingRecHit
TrackingRecHitRefProd refTrkHits
CovarianceMatrix outerStateCovariance() const
outermost trajectory state curvilinear errors
Definition: Track.h:51
unsigned int outerDetId() const
DetId of the detector on which surface the outermost state is located.
Definition: Track.h:59
unsigned int minFound_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
edm::RefProd< TrajectorySeedCollection > refTrajSeeds
edm::ProductID clusterProductB(const TrackingRecHit *hit)
reco::TrackRefProd refTrks
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:230
#define end
Definition: vmac.h:37
std::auto_ptr< reco::TrackExtraCollection > outputTrkExtras
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:63
bool first
Definition: L1TdeRCT.cc:75
edm::AssociationMap< edm::OneToOne< std::vector< Trajectory >, reco::TrackCollection, unsigned short > > TrajTrackAssociationCollection
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:390
TrackListMerger(const edm::ParameterSet &conf)
RefProd< PROD > getRefBeforePut()
Definition: Event.h:128
tuple conf
Definition: dbtoconf.py:185
std::vector< double > indivShareFrac_
bool failedToGet() const
Definition: HandleBase.h:80
reco::TrackExtraRefProd refTrkExtras
std::vector< TrackExtra > TrackExtraCollection
collection of TrackExtra objects
Definition: TrackExtraFwd.h:9
edm::OwnVector< TrackingRecHit > TrackingRecHitCollection
collection of TrackingRecHits
edm::RefToBase< TrajectorySeed > seedRef() const
Definition: Track.h:112
PTrajectoryStateOnDet const & startingState() const
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:46
const math::XYZVector & outerMomentum() const
momentum vector at the outermost hit position
Definition: Track.h:49
bool outerOk() const
return true if the outermost hit is valid
Definition: Track.h:38
std::vector< TkEDGetTokenss > trackProducers_
auto dp
Definition: deltaR.h:24
std::vector< bool > promoteQuality_
edm::Ref< TrackCollection > TrackRef
persistent reference to a Track
Definition: TrackFwd.h:14
CovarianceMatrix innerStateCovariance() const
innermost trajectory state curvilinear errors
Definition: Track.h:53
key_type key() const
Accessor for product key.
Definition: Ref.h:266
range recHits() const
bool isValid() const
void add(const TrackingRecHitRef &r)
add a reference to a RecHit
std::string const & moduleLabel() const
Definition: Provenance.h:60
unsigned int nHits() const
#define likely(x)
Definition: Likely.h:20
#define begin
Definition: vmac.h:30
const math::XYZVector & innerMomentum() const
momentum vector at the innermost hit position
Definition: Track.h:45
TrackingRecHitRef recHit(size_t i) const
Get i-th hit on the track.
Definition: Track.h:67
ProductID id() const
Accessor for product ID.
Definition: RefProd.h:140
const TrackResiduals & residuals() const
Definition: Track.h:117
PropagationDirection seedDirection() const
direction of how the hits were sorted in the original seed
Definition: Track.h:105
tuple cout
Definition: gather_cfg.py:121
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:6
std::auto_ptr< TrackingRecHitCollection > outputTrkHits
ProductID id() const
Accessor for product ID.
Definition: Ref.h:256
DetId geographicalId() const
std::string const & productInstanceName() const
Definition: Provenance.h:62
Provenance getProvenance(BranchID const &theID) const
Definition: Event.cc:76
reco::TrackBase::TrackQuality qualityToSet_
T x() const
Definition: PV3DBase.h:62
virtual LocalPoint localPosition() const =0
unsigned int innerDetId() const
DetId of the detector on which surface the innermost state is located.
Definition: Track.h:61
id_type rawId() const
#define constexpr
void reserve(size_t)
Definition: OwnVector.h:268
edm::RefProd< std::vector< Trajectory > > refTrajs
std::vector< int > hasSelector_
list at
Definition: asciidump.py:428
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:65
std::auto_ptr< TrajTrackAssociationCollection > outputTTAss
void setResiduals(const TrackResiduals &r)
set the residuals
Definition: TrackExtra.h:131