CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
FastTrackMerger.cc
Go to the documentation of this file.
1 #include <memory>
2 
7 
15 
16 
18 
19 #include <vector>
20 #include <map>
21 //
22 
23 //for debug only
24 //#define FAMOS_DEBUG
25 
27 {
28 #ifdef FAMOS_DEBUG
29  std::cout << "FastTrackMerger created" << std::endl;
30 #endif
31 
32  // The main product is a track collection, and all extras
33  produces<reco::TrackCollection>();
34 
35  // The name of the track producers to merge
36  trackProducers = conf.getParameter<std::vector<edm::InputTag> >("TrackProducers");
37 
38  // The name of the track producers to remove
39  std::vector<edm::InputTag> defaultRemove;
41  conf.getUntrackedParameter<std::vector<edm::InputTag> >("RemoveTrackProducers",defaultRemove);
42 
43  // Only the tracks!
44  tracksOnly = conf.getUntrackedParameter<bool>("SaveTracksOnly",false);
45 
46  // optional pT cut
47  double pTMin = conf.getUntrackedParameter<bool>("pTMin",0.);
48  pTMin2 = pTMin*pTMin;
49 
50  // optional nHit cut
51  minHits = conf.getUntrackedParameter<unsigned>("minHits",0);
52 
53  // optional track quality saving
54  promoteQuality = conf.getUntrackedParameter<bool>("promoteTrackQuality",false);
55  qualityStr = conf.getUntrackedParameter<std::string>("newQuality","");
56 
57  // optional trackAlgo (iter0/1/2/3/4)
58  trackAlgo = conf.getUntrackedParameter<unsigned>("trackAlgo",0);
59 
60  //new parameters for Trajectory filtering
61 
62  // The minimum number of hits
63  theMinimumNumberOfHits = conf.getUntrackedParameter<unsigned>("MinNumberOfTrajHits",0);
64 
65  // The maximum number of Lost Hits
66  theMaxLostHits = conf.getUntrackedParameter<unsigned>("MaxLostTrajHits",99);
67 
68  // the max number of consecutive Lost Hits
69  theMaxConsecutiveLostHits = conf.getUntrackedParameter<unsigned>("MaxConsecutiveLostTrajHits",3);
70 
71  //======================
72 
73  if ( !tracksOnly ) {
74  produces<reco::TrackExtraCollection>();
75  produces<TrackingRecHitCollection>();
76  produces<std::vector<Trajectory> >();
77  produces<TrajTrackAssociationCollection>();
78  } else {
79  produces<reco::TrackExtraCollection>();
80  produces<TrackingRecHitCollection>();
81  }
82 
83  // consumes
84  for ( unsigned aProducer=0; aProducer<removeTrackProducers.size(); ++aProducer ) {
85  removeTrackTokens.push_back(consumes<reco::TrackCollection>(removeTrackProducers[aProducer]));
86  }
87  for ( unsigned aProducer=0; aProducer<trackProducers.size(); ++aProducer ) {
88  trackTokens.push_back(consumes<reco::TrackCollection>(trackProducers[aProducer]));
89  trajectoryTokens.push_back(consumes<std::vector<Trajectory> >(trackProducers[aProducer]));
90  assoMapTokens.push_back(consumes<TrajTrackAssociationCollection>(trackProducers[aProducer]));
91  }
92 }
93 
94 
95 // Functions that gets called by framework every event
96 void
98 
99 #ifdef FAMOS_DEBUG
100  std::cout << "################################################################" << std::endl;
101  std::cout << " FastTrackMerger produce init " << std::endl;
102 #endif
103 
104  // The track algorithm (only of merging after iterative tracking)
106  switch(trackAlgo) {
107  case 4: algo = reco::TrackBase::iter0; break;
108  case 5: algo = reco::TrackBase::iter1; break;
109  case 6: algo = reco::TrackBase::iter2; break;
110  case 7: algo = reco::TrackBase::iter3; break;
111  case 8: algo = reco::TrackBase::iter4; break;
112  case 9: algo = reco::TrackBase::iter5; break;
113  case 10: algo = reco::TrackBase::iter6; break;
114  default: algo = reco::TrackBase::undefAlgorithm;
115  }
116 
117  // The produced objects
118  std::auto_ptr<reco::TrackCollection> recoTracks(new reco::TrackCollection);
119  std::auto_ptr<reco::TrackExtraCollection> recoTrackExtras(new reco::TrackExtraCollection);
120  std::auto_ptr<TrackingRecHitCollection> recoHits(new TrackingRecHitCollection);
121  std::auto_ptr<std::vector<Trajectory> > recoTrajectories(new std::vector<Trajectory>);
122  std::auto_ptr<TrajTrackAssociationCollection> recoTrajTrackMap(new TrajTrackAssociationCollection());
123 
124  //prepare the ref to copy the extras
127  std::vector<reco::TrackRef> trackRefs;
128 
129  // No seed -> output an empty track collection
130  if(trackProducers.size() == 0) {
131  e.put(recoTracks);
132  if ( !tracksOnly ) {
133  e.put(recoTrackExtras);
134  e.put(recoHits);
135  e.put(recoTrajectories);
136  e.put(recoTrajTrackMap);
137  } else {
138  e.put(recoTrackExtras);
139  e.put(recoHits);
140  }
141  return;
142  }
143 
144  // The quality to be set
145  reco::TrackBase::TrackQuality qualityToSet;
146  if (qualityStr != "")
148  else
149  qualityToSet = reco::TrackBase::undefQuality;
150 
151  // The input track collection handle
152  edm::Handle<reco::TrackCollection> theTrackCollection;
153 
154  // First, the tracks to be removed
155  std::set<unsigned> removeTracks;
156  for ( unsigned aProducer=0; aProducer<removeTrackProducers.size(); ++aProducer ) {
157  bool isTrackCollection = e.getByToken(removeTrackTokens[aProducer],theTrackCollection);
158  if (!isTrackCollection) continue;
159  reco::TrackCollection::const_iterator aTrack = theTrackCollection->begin();
160  reco::TrackCollection::const_iterator lastTrack = theTrackCollection->end();
161  for ( ; aTrack!=lastTrack; ++aTrack ) {
162  // Get the simtrack Id
163  int recoTrackId = findId(*aTrack);
164  if ( recoTrackId < 0 ) continue;
165  // Remove the track if requested
166  if ( removeTracks.find((unsigned)recoTrackId) != removeTracks.end() ) continue;
167  removeTracks.insert((unsigned)recoTrackId);
168  }
169  }
170 
171  // Then the tracks to be added
172  std::set<unsigned> alreadyAddedTracks;
173 
174  // Loop on the track producers to be merged
175  for ( unsigned aProducer=0; aProducer<trackProducers.size(); ++aProducer ) {
176 
177  bool isTrackCollection = e.getByToken(trackTokens[aProducer],theTrackCollection);
178  if ( ! isTrackCollection ) {
179 #ifdef FAMOS_DEBUG
180  std::cout << "***FastTrackMerger*** Warning! The track collection "
181  << trackProducers[aProducer].encode()
182  << " does not exist." << std::endl;
183 #endif
184  continue;
185  }
186 
187 #ifdef FAMOS_DEBUG
188  std::cout << "***FastTrackMerger*** of track collection "
189  << trackProducers[aProducer].encode()
190  << " with " << theTrackCollection->size()
191  << " tracks to be copied"
192  << std::endl;
193 #endif
194  reco::TrackCollection::const_iterator aTrack = theTrackCollection->begin();
195  reco::TrackCollection::const_iterator lastTrack = theTrackCollection->end();
196 
197  //Copy Tracks and TracksExtra
198  if ( tracksOnly ) {
199 
200  edm:: Handle<reco::TrackExtraCollection > theTrackExtraCollection;
201  // bool isTrackExtraCollection = e.getByLabel(trackProducers[aProducer],theTrackExtraCollection);
202 
203  bool index = 0;
204  for ( ; aTrack!=lastTrack; ++aTrack,++index ) {
205 
206  // Find the track id
207  int recoTrackId = findId(*aTrack);
208  if ( recoTrackId < 0 ) continue;
209 
210  // Ignore tracks to be removed or tracks already copied
211  std::set<unsigned>::iterator iR = removeTracks.find((unsigned)recoTrackId);
212 #ifdef FAMOS_DEBUG
213  if( iR != removeTracks.end() ) std::cout << recoTrackId << "(REMOVED) ";
214 #endif
215  if( iR != removeTracks.end() ) continue;
216 
217  // Ignore tracks already copied
218  std::set<unsigned>::iterator iA = alreadyAddedTracks.find((unsigned)recoTrackId);
219 #ifdef FAMOS_DEBUG
220  if( iA != alreadyAddedTracks.end() ) std::cout << recoTrackId << "(ALREADY ADDED) ";
221 #endif
222  if( iA != alreadyAddedTracks.end() ) continue;
223  //if it is not there then add it!
224  alreadyAddedTracks.insert(recoTrackId);
225 
226 
227 #ifdef FAMOS_DEBUG
228  std::cout << recoTrackId << " ADDED ";
229 #endif
230 
231  // Ignore tracks with too small a pT
232 #ifdef FAMOS_DEBUG
233  if ( aTrack->innerMomentum().Perp2() < pTMin2 )
234  std::cout << "PTMIN CUT APPLIED = " << aTrack->innerMomentum().Perp2() << std::endl;
235 #endif
236 
237  if ( aTrack->innerMomentum().Perp2() < pTMin2 ) continue;
238 
239  // Ignore tracks with too small a pT
240  if ( aTrack->recHitsSize() < minHits ) continue;
241 
242  // A copy of the track + save the transient reference to the track extra reference
243  reco::Track aRecoTrack(*aTrack);
244  recoTracks->push_back(aRecoTrack);
245 
246 
247  recoTrackExtras->push_back(reco::TrackExtra(aRecoTrack.outerPosition(), aRecoTrack.outerMomentum(), aRecoTrack.outerOk(),
248  aRecoTrack.innerPosition(), aRecoTrack.innerMomentum(), aRecoTrack.innerOk(),
249  aRecoTrack.outerStateCovariance(), aRecoTrack.outerDetId(),
250  aRecoTrack.innerStateCovariance(), aRecoTrack.innerDetId(),
251  aRecoTrack.seedDirection(), aRecoTrack.seedRef() ) );
252 
253 
254  recoTracks->back().setExtra(reco::TrackExtraRef(rTrackExtras,recoTrackExtras->size()-1));
255 
256  reco::TrackExtra & tx = recoTrackExtras->back();
257  tx.setResiduals(aRecoTrack.residuals());
258  //TrackingRecHits
259  for( trackingRecHit_iterator hit = aRecoTrack.recHitsBegin(); hit != aRecoTrack.recHitsEnd(); ++ hit ) {
260  recoHits->push_back( (*hit)->clone() );
261  tx.add( TrackingRecHitRef( rHits , recoHits->size() - 1) );
262  }
263 
264  // Save the quality if requested
265  if (promoteQuality) recoTracks->back().setQuality(qualityToSet);
266  if ( trackAlgo ) recoTracks->back().setAlgorithm(algo);
267 
268  }
269 
270 
271  // All extras are to be copied too -> loop on the Trajectory/Track map association
272  } else {
273 
274  edm:: Handle<std::vector<Trajectory> > theTrajectoryCollection;
276 
277  // Count the number of hits and reserve appropriate memory
278  unsigned nRecoHits = 0;
279  for ( ; aTrack!=lastTrack; ++aTrack ) nRecoHits+= aTrack->recHitsSize();
280  recoHits->reserve(nRecoHits); // This is to save some time at push_back.
281 
282  e.getByToken(trajectoryTokens[aProducer],theTrajectoryCollection);
283  e.getByToken(assoMapTokens[aProducer],theAssoMap);
284 
285  // The track collection iterators.
288  anAssociation = theAssoMap->begin();
289  lastAssociation = theAssoMap->end();
290 #ifdef FAMOS_DEBUG
291  std::cout << "List of tracks to be copied " << std::endl;
292 #endif
293  // Build the map of correspondance between reco tracks and sim tracks
294  for ( ; anAssociation != lastAssociation; ++anAssociation ) {
295  edm::Ref<std::vector<Trajectory> > aTrajectoryRef = anAssociation->key;
296  reco::TrackRef aTrackRef = anAssociation->val;
297  // Find the track id
298  int recoTrackId = findId(*aTrackRef);
299  if ( recoTrackId < 0 ) continue;
300 
301  // Ignore tracks to be removed or tracks already copied
302  std::set<unsigned>::iterator iR = removeTracks.find((unsigned)recoTrackId);
303 #ifdef FAMOS_DEBUG
304  if( iR != removeTracks.end() ) std::cout << recoTrackId << "(REMOVED) ";
305 #endif
306  if( iR != removeTracks.end() ) continue;
307 
308  // Ignore tracks already copied
309  std::set<unsigned>::iterator iA = alreadyAddedTracks.find((unsigned)recoTrackId);
310 #ifdef FAMOS_DEBUG
311  if( iA != alreadyAddedTracks.end() ) std::cout << recoTrackId << "(ALREADY ADDED) ";
312 #endif
313  if( iA != alreadyAddedTracks.end() ) continue;
314 
315 #ifdef FAMOS_DEBUG
316  std::cout << recoTrackId << " Newly Added " << std::endl;
317 #endif
318  //if it is not there then add it!
319  alreadyAddedTracks.insert(recoTrackId);
320 
321  // Ignore tracks with too small a pT
322  if ( aTrackRef->innerMomentum().Perp2() < pTMin2 ) continue;
323 
324  // Ignore tracks with too few hits
325  if ( aTrackRef->recHitsSize() < minHits ) continue;
326 
327  //==== add more cuts on the trajectory (emulate the Trajectory Filter)
328 
329 #ifdef FAMOS_DEBUG
330  if(aTrajectoryRef->lostHits() > theMaxLostHits )
331  std::cout << "\tmaxLostHits= " << aTrajectoryRef->lostHits() << "\tCUT =" << theMaxLostHits << std::endl;
332 #endif
333 
334  if((unsigned)aTrajectoryRef->lostHits() > theMaxLostHits ) continue;
335 
336 #ifdef FAMOS_DEBUG
337  if(aTrajectoryRef->foundHits() < theMinimumNumberOfHits )
338  std::cout << "\tMinimumNumberOfHits = " << aTrajectoryRef->foundHits() << "\tCUT = " <<theMinimumNumberOfHits << std::endl;
339 #endif
340 
341  if((unsigned)aTrajectoryRef->foundHits() < theMinimumNumberOfHits ) continue;
342  //calculate the consecutive Lost Hits
343  unsigned consecLostHits = 0;
344  const std::vector<TrajectoryMeasurement> tms = aTrajectoryRef->measurements();
345  for(int itm= tms.size();itm!=0; --itm){
346  if(tms[itm-1].recHit()->isValid())break;
347  else if (Trajectory::lost(*tms[itm-1].recHit())) consecLostHits++;
348  }
349 
350 #ifdef FAMOS_DEBUG
351  if( consecLostHits > theMaxConsecutiveLostHits )
352  std::cout << "\tconsecLostHits = " << consecLostHits << std::endl;
353 #endif
354 
355  if( consecLostHits > theMaxConsecutiveLostHits ) continue;
356 
357 
358  //=============end new filters
359 
360 
361  // A copy of the track
362  reco::Track aRecoTrack(*aTrackRef);
363  recoTracks->push_back(aRecoTrack);
364  // Save the quality if requested
365  if (promoteQuality) recoTracks->back().setQuality(qualityToSet);
366  if ( trackAlgo ) recoTracks->back().setAlgorithm(algo);
367 
368 
369  // A copy of the hits
370  unsigned nh = aRecoTrack.recHitsSize();
371  for ( unsigned ih=0; ih<nh; ++ih ) {
372  TrackingRecHit *hit = aRecoTrack.recHit(ih)->clone();
373  recoHits->push_back(hit);
374  }
375 
376  // A copy of the trajectories
377  recoTrajectories->push_back(*aTrajectoryRef);
378 
379  }
380 #ifdef FAMOS_DEBUG
381  std::cout << std::endl;
382 #endif
383  }
384 }
385 
386  // Save the track candidates in the event
387 #ifdef FAMOS_DEBUG
388  std::cout << "Saving "
389  << recoTracks->size() << " reco::Tracks " << std::endl;
390 
391 #endif
392 
393  if ( tracksOnly ) {
394  // Save only the tracks (with transient reference to track extras)
395  e.put(recoTracks);
396  e.put(recoTrackExtras);
397  e.put(recoHits);
398 
399 
400  } else {
401  // Save the tracking recHits
402  edm::OrphanHandle <TrackingRecHitCollection> theRecoHits = e.put(recoHits );
403 
404  // Create the track extras and add the references to the rechits
405  unsigned hits=0;
406  unsigned nTracks = recoTracks->size();
407  recoTrackExtras->reserve(nTracks); // To save some time at push_back
408  for ( unsigned index = 0; index < nTracks; ++index ) {
409 
410  //reco::TrackExtra aTrackExtra;
411  reco::Track& aTrack = recoTracks->at(index);
412 
413  // std::cout << "initial TrackAlgo = " << trackAlgo << "\tNAME " << aTrack.algoName() << "\tnumber = " << aTrack.algo() << std::endl;
414 
415  reco::TrackExtra aTrackExtra(aTrack.outerPosition(),
416  aTrack.outerMomentum(),
417  aTrack.outerOk(),
418  aTrack.innerPosition(),
419  aTrack.innerMomentum(),
420  aTrack.innerOk(),
421  aTrack.outerStateCovariance(),
422  aTrack.outerDetId(),
423  aTrack.innerStateCovariance(),
424  aTrack.innerDetId(),
425  aTrack.seedDirection(),
426  aTrack.seedRef());
427 
428  unsigned nHits = aTrack.recHitsSize();
429  for ( unsigned int ih=0; ih<nHits; ++ih) {
430  aTrackExtra.add(TrackingRecHitRef(theRecoHits,hits++));
431  }
432  recoTrackExtras->push_back(aTrackExtra);
433  }
434 
435  // Save the track extras
436  edm::OrphanHandle<reco::TrackExtraCollection> theRecoTrackExtras = e.put(recoTrackExtras);
437 
438  // Add the reference to the track extra in the tracks
439  for ( unsigned index = 0; index<nTracks; ++index ) {
440  const reco::TrackExtraRef theTrackExtraRef(theRecoTrackExtras,index);
441  (recoTracks->at(index)).setExtra(theTrackExtraRef);
442  }
443 
444  // Save the tracks
445  edm::OrphanHandle<reco::TrackCollection> theRecoTracks = e.put(recoTracks);
446 
447  // Save the trajectories
448  edm::OrphanHandle<std::vector<Trajectory> > theRecoTrajectories = e.put(recoTrajectories);
449 
450  // Create and set the trajectory/track association map
451  for ( unsigned index = 0; index<nTracks; ++index ) {
452  edm::Ref<std::vector<Trajectory> > trajRef( theRecoTrajectories, index );
453  edm::Ref<reco::TrackCollection> tkRef( theRecoTracks, index );
454  recoTrajTrackMap->insert(trajRef,tkRef);
455  }
456 
457  // Save the association map.
458  e.put(recoTrajTrackMap);
459 
460  }
461 
462 }
463 
464 int
465 FastTrackMerger::findId(const reco::Track& aTrack) const {
466  int trackId = -1;
467  trackingRecHit_iterator aHit = aTrack.recHitsBegin();
468  trackingRecHit_iterator lastHit = aTrack.recHitsEnd();
469  for ( ; aHit!=lastHit; ++aHit ) {
470  if ( !aHit->get()->isValid() ) continue;
471  // const SiTrackerGSRecHit2D * rechit = (const SiTrackerGSRecHit2D*) (aHit->get());
472  const SiTrackerGSMatchedRecHit2D * rechit = (const SiTrackerGSMatchedRecHit2D*) (aHit->get());
473  trackId = rechit->simtrackId();
474  break;
475  }
476  return trackId;
477 }
478 
479 
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
unsigned theMinimumNumberOfHits
std::vector< edm::EDGetTokenT< reco::TrackCollection > > removeTrackTokens
std::vector< edm::InputTag > removeTrackProducers
int findId(const reco::Track &aTrack) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
FastTrackMerger(const edm::ParameterSet &conf)
TrackQuality
track quality
Definition: TrackBase.h:113
size_t recHitsSize() const
Get number of RecHits. (Warning, this includes invalid hits, which are not physical hits)...
Definition: Track.h:119
std::vector< edm::EDGetTokenT< std::vector< Trajectory > > > trajectoryTokens
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:13
bool innerOk() const
return true if the innermost hit is valid
Definition: Track.h:50
key_type key() const
Accessor for product key.
Definition: Ref.h:266
const math::XYZPoint & outerPosition() const
position of the outermost hit
Definition: Track.h:65
TrackAlgorithm
track algorithm
Definition: TrackBase.h:97
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
const math::XYZPoint & innerPosition() const
position of the innermost hit
Definition: Track.h:55
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
edm::Ref< TrackingRecHitCollection > TrackingRecHitRef
persistent reference to a TrackingRecHit
CovarianceMatrix outerStateCovariance() const
outermost trajectory state curvilinear errors
Definition: Track.h:75
unsigned int outerDetId() const
DetId of the detector on which surface the outermost state is located.
Definition: Track.h:94
std::vector< edm::EDGetTokenT< reco::TrackCollection > > trackTokens
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:104
edm::AssociationMap< edm::OneToOne< std::vector< Trajectory >, reco::TrackCollection, unsigned short > > TrajTrackAssociationCollection
RefProd< PROD > getRefBeforePut()
Definition: Event.h:133
tuple conf
Definition: dbtoconf.py:185
std::string qualityStr
unsigned theMaxConsecutiveLostHits
std::vector< TrackExtra > TrackExtraCollection
collection of TrackExtra objects
Definition: TrackExtraFwd.h:11
edm::RefToBase< TrajectorySeed > seedRef() const
Definition: Track.h:213
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:102
const math::XYZVector & outerMomentum() const
momentum vector at the outermost hit position
Definition: Track.h:70
bool outerOk() const
return true if the outermost hit is valid
Definition: Track.h:45
virtual void produce(edm::Event &e, const edm::EventSetup &es)
CovarianceMatrix innerStateCovariance() const
innermost trajectory state curvilinear errors
Definition: Track.h:80
void add(const TrackingRecHitRef &r)
add a reference to a RecHit
std::vector< edm::InputTag > trackProducers
std::vector< edm::EDGetTokenT< TrajTrackAssociationCollection > > assoMapTokens
const math::XYZVector & innerMomentum() const
momentum vector at the innermost hit position
Definition: Track.h:60
TrackingRecHitRef recHit(size_t i) const
Get i-th hit on the track.
Definition: Track.h:114
const TrackResiduals & residuals() const
Definition: Track.h:220
PropagationDirection seedDirection() const
direction of how the hits were sorted in the original seed
Definition: Track.h:204
unsigned theMaxLostHits
tuple cout
Definition: gather_cfg.py:121
unsigned int innerDetId() const
DetId of the detector on which surface the innermost state is located.
Definition: Track.h:99
static bool lost(const TrackingRecHit &hit)
Definition: Trajectory.cc:128
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:109
void setResiduals(const TrackResiduals &r)
set the residuals
Definition: TrackExtra.h:178