CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
FastTrackMerger Class Reference

#include <FastTrackMerger.h>

Inheritance diagram for FastTrackMerger:
edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

Public Member Functions

 FastTrackMerger (const edm::ParameterSet &conf)
 
virtual void produce (edm::Event &e, const edm::EventSetup &es)
 
virtual ~FastTrackMerger ()
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
boost::function< void(const
BranchDescription &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 

Private Member Functions

int findId (const reco::Track &aTrack) const
 

Private Attributes

unsigned minHits
 
bool promoteQuality
 
double pTMin2
 
std::string qualityStr
 
std::vector< edm::InputTagremoveTrackProducers
 
unsigned theMaxConsecutiveLostHits
 
unsigned theMaxLostHits
 
unsigned theMinimumNumberOfHits
 
unsigned trackAlgo
 
std::vector< edm::InputTagtrackProducers
 
bool tracksOnly
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
typedef WorkerT< EDProducerWorkerType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDProducer
CurrentProcessingContext const * currentContext () const
 
- Protected Member Functions inherited from edm::ProducerBase
template<class TProducer , class TMethod >
void callWhenNewProductsRegistered (TProducer *iProd, TMethod iMethod)
 

Detailed Description

Definition at line 20 of file FastTrackMerger.h.

Constructor & Destructor Documentation

FastTrackMerger::FastTrackMerger ( const edm::ParameterSet conf)
explicit

Definition at line 29 of file FastTrackMerger.cc.

References gather_cfg::cout, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), minHits, promoteQuality, pTMin2, qualityStr, removeTrackProducers, theMaxConsecutiveLostHits, theMaxLostHits, theMinimumNumberOfHits, trackAlgo, trackProducers, and tracksOnly.

30 {
31 #ifdef FAMOS_DEBUG
32  std::cout << "FastTrackMerger created" << std::endl;
33 #endif
34 
35  // The main product is a track collection, and all extras
36  produces<reco::TrackCollection>();
37 
38  // The name of the track producers to merge
39  trackProducers = conf.getParameter<std::vector<edm::InputTag> >("TrackProducers");
40 
41  // The name of the track producers to remove
42  std::vector<edm::InputTag> defaultRemove;
44  conf.getUntrackedParameter<std::vector<edm::InputTag> >("RemoveTrackProducers",defaultRemove);
45 
46  // Only the tracks!
47  tracksOnly = conf.getUntrackedParameter<bool>("SaveTracksOnly",false);
48 
49  // optional pT cut
50  double pTMin = conf.getUntrackedParameter<bool>("pTMin",0.);
51  pTMin2 = pTMin*pTMin;
52 
53  // optional nHit cut
54  minHits = conf.getUntrackedParameter<unsigned>("minHits",0);
55 
56  // optional track quality saving
57  promoteQuality = conf.getUntrackedParameter<bool>("promoteTrackQuality",false);
58  qualityStr = conf.getUntrackedParameter<std::string>("newQuality","");
59 
60  // optional trackAlgo (iter0/1/2/3/4)
61  trackAlgo = conf.getUntrackedParameter<unsigned>("trackAlgo",0);
62 
63  //new parameters for Trajectory filtering
64 
65  // The minimum number of hits
66  theMinimumNumberOfHits = conf.getUntrackedParameter<unsigned>("MinNumberOfTrajHits",0);
67 
68  // The maximum number of Lost Hits
69  theMaxLostHits = conf.getUntrackedParameter<unsigned>("MaxLostTrajHits",99);
70 
71  // the max number of consecutive Lost Hits
72  theMaxConsecutiveLostHits = conf.getUntrackedParameter<unsigned>("MaxConsecutiveLostTrajHits",3);
73 
74  //======================
75 
76  if ( !tracksOnly ) {
77  produces<reco::TrackExtraCollection>();
78  produces<TrackingRecHitCollection>();
79  produces<std::vector<Trajectory> >();
80  produces<TrajTrackAssociationCollection>();
81  } else {
82  produces<reco::TrackExtraCollection>();
83  produces<TrackingRecHitCollection>();
84  }
85 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
unsigned theMinimumNumberOfHits
std::vector< edm::InputTag > removeTrackProducers
std::string qualityStr
unsigned theMaxConsecutiveLostHits
std::vector< edm::InputTag > trackProducers
unsigned theMaxLostHits
tuple cout
Definition: gather_cfg.py:41
virtual FastTrackMerger::~FastTrackMerger ( )
inlinevirtual

Definition at line 26 of file FastTrackMerger.h.

26 {}

Member Function Documentation

int FastTrackMerger::findId ( const reco::Track aTrack) const
private

Definition at line 457 of file FastTrackMerger.cc.

References reco::Track::recHitsBegin(), reco::Track::recHitsEnd(), and SiTrackerGSMatchedRecHit2D::simtrackId().

Referenced by produce().

457  {
458  int trackId = -1;
459  trackingRecHit_iterator aHit = aTrack.recHitsBegin();
460  trackingRecHit_iterator lastHit = aTrack.recHitsEnd();
461  for ( ; aHit!=lastHit; ++aHit ) {
462  if ( !aHit->get()->isValid() ) continue;
463  // const SiTrackerGSRecHit2D * rechit = (const SiTrackerGSRecHit2D*) (aHit->get());
464  const SiTrackerGSMatchedRecHit2D * rechit = (const SiTrackerGSMatchedRecHit2D*) (aHit->get());
465  trackId = rechit->simtrackId();
466  break;
467  }
468  return trackId;
469 }
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:63
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:65
void FastTrackMerger::produce ( edm::Event e,
const edm::EventSetup es 
)
virtual

Implements edm::EDProducer.

Definition at line 90 of file FastTrackMerger.cc.

References reco::TrackExtraBase::add(), gather_cfg::cout, findId(), edm::Event::getByLabel(), edm::Event::getRefBeforePut(), getHLTprescales::index, reco::Track::innerDetId(), reco::Track::innerMomentum(), reco::Track::innerOk(), reco::Track::innerPosition(), reco::Track::innerStateCovariance(), reco::TrackBase::iter0, reco::TrackBase::iter1, reco::TrackBase::iter2, reco::TrackBase::iter3, reco::TrackBase::iter4, reco::TrackBase::iter5, edm::helpers::KeyVal< K, V >::key, Trajectory::lost(), minHits, reco::Track::outerDetId(), reco::Track::outerMomentum(), reco::Track::outerOk(), reco::Track::outerPosition(), reco::Track::outerStateCovariance(), promoteQuality, pTMin2, edm::Event::put(), reco::TrackBase::qualityByName(), qualityStr, reco::Track::recHit(), reco::Track::recHitsBegin(), reco::Track::recHitsEnd(), reco::Track::recHitsSize(), removeTrackProducers, reco::Track::residuals(), reco::Track::seedDirection(), reco::Track::seedRef(), reco::TrackExtra::setResiduals(), theMaxConsecutiveLostHits, theMaxLostHits, theMinimumNumberOfHits, trackAlgo, trackProducers, tracksOnly, reco::TrackBase::undefAlgorithm, reco::TrackBase::undefQuality, and edm::helpers::KeyVal< K, V >::val.

Referenced by python.JSONExport.JsonExport::export(), and python.HTMLExport.HTMLExport::export().

90  {
91 
92 #ifdef FAMOS_DEBUG
93  std::cout << "################################################################" << std::endl;
94  std::cout << " FastTrackMerger produce init " << std::endl;
95 #endif
96 
97  // The track algorithm (only of merging after iterative tracking)
99  switch(trackAlgo) {
100  case 4: algo = reco::TrackBase::iter0; break;
101  case 5: algo = reco::TrackBase::iter1; break;
102  case 6: algo = reco::TrackBase::iter2; break;
103  case 7: algo = reco::TrackBase::iter3; break;
104  case 8: algo = reco::TrackBase::iter4; break;
105  case 9: algo = reco::TrackBase::iter5; break;
106  default: algo = reco::TrackBase::undefAlgorithm;
107  }
108 
109  // The produced objects
110  std::auto_ptr<reco::TrackCollection> recoTracks(new reco::TrackCollection);
111  std::auto_ptr<reco::TrackExtraCollection> recoTrackExtras(new reco::TrackExtraCollection);
112  std::auto_ptr<TrackingRecHitCollection> recoHits(new TrackingRecHitCollection);
113  std::auto_ptr<std::vector<Trajectory> > recoTrajectories(new std::vector<Trajectory>);
114  std::auto_ptr<TrajTrackAssociationCollection> recoTrajTrackMap(new TrajTrackAssociationCollection());
115 
116  //prepare the ref to copy the extras
119  std::vector<reco::TrackRef> trackRefs;
120 
121  // No seed -> output an empty track collection
122  if(trackProducers.size() == 0) {
123  e.put(recoTracks);
124  if ( !tracksOnly ) {
125  e.put(recoTrackExtras);
126  e.put(recoHits);
127  e.put(recoTrajectories);
128  e.put(recoTrajTrackMap);
129  } else {
130  e.put(recoTrackExtras);
131  e.put(recoHits);
132  }
133  return;
134  }
135 
136  // The quality to be set
137  reco::TrackBase::TrackQuality qualityToSet;
138  if (qualityStr != "")
140  else
141  qualityToSet = reco::TrackBase::undefQuality;
142 
143  // The input track collection handle
144  edm::Handle<reco::TrackCollection> theTrackCollection;
145 
146  // First, the tracks to be removed
147  std::set<unsigned> removeTracks;
148  for ( unsigned aProducer=0; aProducer<removeTrackProducers.size(); ++aProducer ) {
149  bool isTrackCollection = e.getByLabel(removeTrackProducers[aProducer],theTrackCollection);
150  if (!isTrackCollection) continue;
151  reco::TrackCollection::const_iterator aTrack = theTrackCollection->begin();
152  reco::TrackCollection::const_iterator lastTrack = theTrackCollection->end();
153  for ( ; aTrack!=lastTrack; ++aTrack ) {
154  // Get the simtrack Id
155  int recoTrackId = findId(*aTrack);
156  if ( recoTrackId < 0 ) continue;
157  // Remove the track if requested
158  if ( removeTracks.find((unsigned)recoTrackId) != removeTracks.end() ) continue;
159  removeTracks.insert((unsigned)recoTrackId);
160  }
161  }
162 
163  // Then the tracks to be added
164  std::set<unsigned> alreadyAddedTracks;
165 
166  // Loop on the track producers to be merged
167  for ( unsigned aProducer=0; aProducer<trackProducers.size(); ++aProducer ) {
168 
169  bool isTrackCollection = e.getByLabel(trackProducers[aProducer],theTrackCollection);
170  if ( ! isTrackCollection ) {
171 #ifdef FAMOS_DEBUG
172  std::cout << "***FastTrackMerger*** Warning! The track collection "
173  << trackProducers[aProducer].encode()
174  << " does not exist." << std::endl;
175 #endif
176  continue;
177  }
178 
179 #ifdef FAMOS_DEBUG
180  std::cout << "***FastTrackMerger*** of track collection "
181  << trackProducers[aProducer].encode()
182  << " with " << theTrackCollection->size()
183  << " tracks to be copied"
184  << std::endl;
185 #endif
186  reco::TrackCollection::const_iterator aTrack = theTrackCollection->begin();
187  reco::TrackCollection::const_iterator lastTrack = theTrackCollection->end();
188 
189  //Copy Tracks and TracksExtra
190  if ( tracksOnly ) {
191 
192  edm:: Handle<reco::TrackExtraCollection > theTrackExtraCollection;
193  // bool isTrackExtraCollection = e.getByLabel(trackProducers[aProducer],theTrackExtraCollection);
194 
195  bool index = 0;
196  for ( ; aTrack!=lastTrack; ++aTrack,++index ) {
197 
198  // Find the track id
199  int recoTrackId = findId(*aTrack);
200  if ( recoTrackId < 0 ) continue;
201 
202  // Ignore tracks to be removed or tracks already copied
203  std::set<unsigned>::iterator iR = removeTracks.find((unsigned)recoTrackId);
204 #ifdef FAMOS_DEBUG
205  if( iR != removeTracks.end() ) std::cout << recoTrackId << "(REMOVED) ";
206 #endif
207  if( iR != removeTracks.end() ) continue;
208 
209  // Ignore tracks already copied
210  std::set<unsigned>::iterator iA = alreadyAddedTracks.find((unsigned)recoTrackId);
211 #ifdef FAMOS_DEBUG
212  if( iA != alreadyAddedTracks.end() ) std::cout << recoTrackId << "(ALREADY ADDED) ";
213 #endif
214  if( iA != alreadyAddedTracks.end() ) continue;
215  //if it is not there then add it!
216  alreadyAddedTracks.insert(recoTrackId);
217 
218 
219 #ifdef FAMOS_DEBUG
220  std::cout << recoTrackId << " ADDED ";
221 #endif
222 
223  // Ignore tracks with too small a pT
224 #ifdef FAMOS_DEBUG
225  if ( aTrack->innerMomentum().Perp2() < pTMin2 )
226  std::cout << "PTMIN CUT APPLIED = " << aTrack->innerMomentum().Perp2() << std::endl;
227 #endif
228 
229  if ( aTrack->innerMomentum().Perp2() < pTMin2 ) continue;
230 
231  // Ignore tracks with too small a pT
232  if ( aTrack->recHitsSize() < minHits ) continue;
233 
234  // A copy of the track + save the transient reference to the track extra reference
235  reco::Track aRecoTrack(*aTrack);
236  recoTracks->push_back(aRecoTrack);
237 
238 
239  recoTrackExtras->push_back(reco::TrackExtra(aRecoTrack.outerPosition(), aRecoTrack.outerMomentum(), aRecoTrack.outerOk(),
240  aRecoTrack.innerPosition(), aRecoTrack.innerMomentum(), aRecoTrack.innerOk(),
241  aRecoTrack.outerStateCovariance(), aRecoTrack.outerDetId(),
242  aRecoTrack.innerStateCovariance(), aRecoTrack.innerDetId(),
243  aRecoTrack.seedDirection(), aRecoTrack.seedRef() ) );
244 
245 
246  recoTracks->back().setExtra(reco::TrackExtraRef(rTrackExtras,recoTrackExtras->size()-1));
247 
248  reco::TrackExtra & tx = recoTrackExtras->back();
249  tx.setResiduals(aRecoTrack.residuals());
250  //TrackingRecHits
251  for( trackingRecHit_iterator hit = aRecoTrack.recHitsBegin(); hit != aRecoTrack.recHitsEnd(); ++ hit ) {
252  recoHits->push_back( (*hit)->clone() );
253  tx.add( TrackingRecHitRef( rHits , recoHits->size() - 1) );
254  }
255 
256  // Save the quality if requested
257  if (promoteQuality) recoTracks->back().setQuality(qualityToSet);
258  if ( trackAlgo ) recoTracks->back().setAlgorithm(algo);
259 
260  }
261 
262 
263  // All extras are to be copied too -> loop on the Trajectory/Track map association
264  } else {
265 
266  edm:: Handle<std::vector<Trajectory> > theTrajectoryCollection;
268 
269  // Count the number of hits and reserve appropriate memory
270  unsigned nRecoHits = 0;
271  for ( ; aTrack!=lastTrack; ++aTrack ) nRecoHits+= aTrack->recHitsSize();
272  recoHits->reserve(nRecoHits); // This is to save some time at push_back.
273 
274  e.getByLabel(trackProducers[aProducer],theTrajectoryCollection);
275  e.getByLabel(trackProducers[aProducer],theAssoMap);
276 
277  // The track collection iterators.
280  anAssociation = theAssoMap->begin();
281  lastAssociation = theAssoMap->end();
282 #ifdef FAMOS_DEBUG
283  std::cout << "List of tracks to be copied " << std::endl;
284 #endif
285  // Build the map of correspondance between reco tracks and sim tracks
286  for ( ; anAssociation != lastAssociation; ++anAssociation ) {
287  edm::Ref<std::vector<Trajectory> > aTrajectoryRef = anAssociation->key;
288  reco::TrackRef aTrackRef = anAssociation->val;
289  // Find the track id
290  int recoTrackId = findId(*aTrackRef);
291  if ( recoTrackId < 0 ) continue;
292 
293  // Ignore tracks to be removed or tracks already copied
294  std::set<unsigned>::iterator iR = removeTracks.find((unsigned)recoTrackId);
295 #ifdef FAMOS_DEBUG
296  if( iR != removeTracks.end() ) std::cout << recoTrackId << "(REMOVED) ";
297 #endif
298  if( iR != removeTracks.end() ) continue;
299 
300  // Ignore tracks already copied
301  std::set<unsigned>::iterator iA = alreadyAddedTracks.find((unsigned)recoTrackId);
302 #ifdef FAMOS_DEBUG
303  if( iA != alreadyAddedTracks.end() ) std::cout << recoTrackId << "(ALREADY ADDED) ";
304 #endif
305  if( iA != alreadyAddedTracks.end() ) continue;
306 
307 #ifdef FAMOS_DEBUG
308  std::cout << recoTrackId << " Newly Added " << std::endl;
309 #endif
310  //if it is not there then add it!
311  alreadyAddedTracks.insert(recoTrackId);
312 
313  // Ignore tracks with too small a pT
314  if ( aTrackRef->innerMomentum().Perp2() < pTMin2 ) continue;
315 
316  // Ignore tracks with too few hits
317  if ( aTrackRef->recHitsSize() < minHits ) continue;
318 
319  //==== add more cuts on the trajectory (emulate the Trajectory Filter)
320 
321 #ifdef FAMOS_DEBUG
322  if(aTrajectoryRef->lostHits() > theMaxLostHits )
323  std::cout << "\tmaxLostHits= " << aTrajectoryRef->lostHits() << "\tCUT =" << theMaxLostHits << std::endl;
324 #endif
325 
326  if((unsigned)aTrajectoryRef->lostHits() > theMaxLostHits ) continue;
327 
328 #ifdef FAMOS_DEBUG
329  if(aTrajectoryRef->foundHits() < theMinimumNumberOfHits )
330  std::cout << "\tMinimumNumberOfHits = " << aTrajectoryRef->foundHits() << "\tCUT = " <<theMinimumNumberOfHits << std::endl;
331 #endif
332 
333  if((unsigned)aTrajectoryRef->foundHits() < theMinimumNumberOfHits ) continue;
334  //calculate the consecutive Lost Hits
335  unsigned consecLostHits = 0;
336  const std::vector<TrajectoryMeasurement> tms = aTrajectoryRef->measurements();
337  for(int itm= tms.size();itm!=0; --itm){
338  if(tms[itm-1].recHit()->isValid())break;
339  else if (Trajectory::lost(*tms[itm-1].recHit())) consecLostHits++;
340  }
341 
342 #ifdef FAMOS_DEBUG
343  if( consecLostHits > theMaxConsecutiveLostHits )
344  std::cout << "\tconsecLostHits = " << consecLostHits << std::endl;
345 #endif
346 
347  if( consecLostHits > theMaxConsecutiveLostHits ) continue;
348 
349 
350  //=============end new filters
351 
352 
353  // A copy of the track
354  reco::Track aRecoTrack(*aTrackRef);
355  recoTracks->push_back(aRecoTrack);
356  // Save the quality if requested
357  if (promoteQuality) recoTracks->back().setQuality(qualityToSet);
358  if ( trackAlgo ) recoTracks->back().setAlgorithm(algo);
359 
360 
361  // A copy of the hits
362  unsigned nh = aRecoTrack.recHitsSize();
363  for ( unsigned ih=0; ih<nh; ++ih ) {
364  TrackingRecHit *hit = aRecoTrack.recHit(ih)->clone();
365  recoHits->push_back(hit);
366  }
367 
368  // A copy of the trajectories
369  recoTrajectories->push_back(*aTrajectoryRef);
370 
371  }
372 #ifdef FAMOS_DEBUG
373  std::cout << std::endl;
374 #endif
375  }
376 }
377 
378  // Save the track candidates in the event
379 #ifdef FAMOS_DEBUG
380  std::cout << "Saving "
381  << recoTracks->size() << " reco::Tracks " << std::endl;
382 
383 #endif
384 
385  if ( tracksOnly ) {
386  // Save only the tracks (with transient reference to track extras)
387  e.put(recoTracks);
388  e.put(recoTrackExtras);
389  e.put(recoHits);
390 
391 
392  } else {
393  // Save the tracking recHits
394  edm::OrphanHandle <TrackingRecHitCollection> theRecoHits = e.put(recoHits );
395 
396  // Create the track extras and add the references to the rechits
397  unsigned hits=0;
398  unsigned nTracks = recoTracks->size();
399  recoTrackExtras->reserve(nTracks); // To save some time at push_back
400  for ( unsigned index = 0; index < nTracks; ++index ) {
401 
402  //reco::TrackExtra aTrackExtra;
403  reco::Track& aTrack = recoTracks->at(index);
404 
405  // std::cout << "initial TrackAlgo = " << trackAlgo << "\tNAME " << aTrack.algoName() << "\tnumber = " << aTrack.algo() << std::endl;
406 
407  reco::TrackExtra aTrackExtra(aTrack.outerPosition(),
408  aTrack.outerMomentum(),
409  aTrack.outerOk(),
410  aTrack.innerPosition(),
411  aTrack.innerMomentum(),
412  aTrack.innerOk(),
413  aTrack.outerStateCovariance(),
414  aTrack.outerDetId(),
415  aTrack.innerStateCovariance(),
416  aTrack.innerDetId(),
417  aTrack.seedDirection(),
418  aTrack.seedRef());
419 
420  unsigned nHits = aTrack.recHitsSize();
421  for ( unsigned int ih=0; ih<nHits; ++ih) {
422  aTrackExtra.add(TrackingRecHitRef(theRecoHits,hits++));
423  }
424  recoTrackExtras->push_back(aTrackExtra);
425  }
426 
427  // Save the track extras
428  edm::OrphanHandle<reco::TrackExtraCollection> theRecoTrackExtras = e.put(recoTrackExtras);
429 
430  // Add the reference to the track extra in the tracks
431  for ( unsigned index = 0; index<nTracks; ++index ) {
432  const reco::TrackExtraRef theTrackExtraRef(theRecoTrackExtras,index);
433  (recoTracks->at(index)).setExtra(theTrackExtraRef);
434  }
435 
436  // Save the tracks
437  edm::OrphanHandle<reco::TrackCollection> theRecoTracks = e.put(recoTracks);
438 
439  // Save the trajectories
440  edm::OrphanHandle<std::vector<Trajectory> > theRecoTrajectories = e.put(recoTrajectories);
441 
442  // Create and set the trajectory/track association map
443  for ( unsigned index = 0; index<nTracks; ++index ) {
444  edm::Ref<std::vector<Trajectory> > trajRef( theRecoTrajectories, index );
445  edm::Ref<reco::TrackCollection> tkRef( theRecoTracks, index );
446  recoTrajTrackMap->insert(trajRef,tkRef);
447  }
448 
449  // Save the association map.
450  e.put(recoTrajTrackMap);
451 
452  }
453 
454 }
unsigned theMinimumNumberOfHits
static bool lost(const TransientTrackingRecHit &hit)
Definition: Trajectory.cc:205
std::vector< edm::InputTag > removeTrackProducers
int findId(const reco::Track &aTrack) const
TrackQuality
track quality
Definition: TrackBase.h:95
size_t recHitsSize() const
Get number of RecHits. (Warning, this includes invalid hits, which are not physical hits)...
Definition: Track.h:69
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:10
bool innerOk() const
return true if the innermost hit is valid
Definition: Track.h:40
const math::XYZPoint & outerPosition() const
position of the outermost hit
Definition: Track.h:47
TrackAlgorithm
track algorithm
Definition: TrackBase.h:82
const math::XYZPoint & innerPosition() const
position of the innermost hit
Definition: Track.h:42
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:84
edm::Ref< TrackingRecHitCollection > TrackingRecHitRef
persistent reference to a TrackingRecHit
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
edm::AssociationMap< edm::OneToOne< std::vector< Trajectory >, reco::TrackCollection, unsigned short > > TrajTrackAssociationCollection
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:355
RefProd< PROD > getRefBeforePut()
Definition: Event.h:96
std::string qualityStr
unsigned theMaxConsecutiveLostHits
std::vector< TrackExtra > TrackExtraCollection
collection of TrackExtra objects
Definition: TrackExtraFwd.h:9
edm::RefToBase< TrajectorySeed > seedRef() const
Definition: Track.h:112
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
CovarianceMatrix innerStateCovariance() const
innermost trajectory state curvilinear errors
Definition: Track.h:53
void add(const TrackingRecHitRef &r)
add a reference to a RecHit
std::vector< edm::InputTag > trackProducers
const math::XYZVector & innerMomentum() const
momentum vector at the innermost hit position
Definition: Track.h:45
PropagationDirection seedDirection() const
direction of how the hits were sorted in the original seed
Definition: Track.h:105
unsigned theMaxLostHits
tuple cout
Definition: gather_cfg.py:41
unsigned int innerDetId() const
DetId of the detector on which surface the innermost state is located.
Definition: Track.h:61
void setResiduals(const TrackResiduals &r)
set the residuals
Definition: TrackExtra.h:131

Member Data Documentation

unsigned FastTrackMerger::minHits
private

Definition at line 41 of file FastTrackMerger.h.

Referenced by FastTrackMerger(), and produce().

bool FastTrackMerger::promoteQuality
private

Definition at line 39 of file FastTrackMerger.h.

Referenced by FastTrackMerger(), and produce().

double FastTrackMerger::pTMin2
private

Definition at line 40 of file FastTrackMerger.h.

Referenced by FastTrackMerger(), and produce().

std::string FastTrackMerger::qualityStr
private

Definition at line 43 of file FastTrackMerger.h.

Referenced by FastTrackMerger(), and produce().

std::vector<edm::InputTag> FastTrackMerger::removeTrackProducers
private

Definition at line 37 of file FastTrackMerger.h.

Referenced by FastTrackMerger(), and produce().

unsigned FastTrackMerger::theMaxConsecutiveLostHits
private

Definition at line 46 of file FastTrackMerger.h.

Referenced by FastTrackMerger(), and produce().

unsigned FastTrackMerger::theMaxLostHits
private

Definition at line 45 of file FastTrackMerger.h.

Referenced by FastTrackMerger(), and produce().

unsigned FastTrackMerger::theMinimumNumberOfHits
private

Definition at line 44 of file FastTrackMerger.h.

Referenced by FastTrackMerger(), and produce().

unsigned FastTrackMerger::trackAlgo
private

Definition at line 42 of file FastTrackMerger.h.

Referenced by FastTrackMerger(), and produce().

std::vector<edm::InputTag> FastTrackMerger::trackProducers
private

Definition at line 36 of file FastTrackMerger.h.

Referenced by FastTrackMerger(), and produce().

bool FastTrackMerger::tracksOnly
private

Definition at line 38 of file FastTrackMerger.h.

Referenced by FastTrackMerger(), and produce().