CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrackCandidateProducer.cc
Go to the documentation of this file.
1 #include <memory>
2 
7 
17 
20 
23 
25 //#include "FastSimulation/Tracking/interface/TrackerRecHitSplit.h"
26 
28 
29 #include <vector>
30 #include <map>
31 
36 
42 
43 //Propagator withMaterial
45 
48 
49 //
50 
51 //for debug only
52 //#define FAMOS_DEBUG
53 
55 {
56 #ifdef FAMOS_DEBUG
57  std::cout << "TrackCandidateProducer created" << std::endl;
58 #endif
59 
60  // The main product is a track candidate collection.
61  produces<TrackCandidateCollection>();
62 
63  // These products contain tracks already reconstructed at this level
64  // (No need to reconstruct them twice!)
65  produces<reco::TrackCollection>();
66  produces<TrackingRecHitCollection>();
67  produces<reco::TrackExtraCollection>();
68  produces<std::vector<Trajectory> >();
69  produces<TrajTrackAssociationCollection>();
70 
71  // The name of the seed producer
72  seedProducer = conf.getParameter<edm::InputTag>("SeedProducer");
73 
74  // The name of the recHit producer
75  hitProducer = conf.getParameter<edm::InputTag>("HitProducer");
76 
77  // The name of the track producer (tracks already produced need not be produced again!)
78  // trackProducer = conf.getParameter<edm::InputTag>("TrackProducer");
79  trackProducers = conf.getParameter<std::vector<edm::InputTag> >("TrackProducers");
80 
81  // Copy (or not) the tracks already produced in a new collection
82  keepFittedTracks = conf.getParameter<bool>("KeepFittedTracks");
83 
84  // The minimum number of crossed layers
85  minNumberOfCrossedLayers = conf.getParameter<unsigned int>("MinNumberOfCrossedLayers");
86 
87  // The maximum number of crossed layers
88  maxNumberOfCrossedLayers = conf.getParameter<unsigned int>("MaxNumberOfCrossedLayers");
89 
90  // Reject overlapping hits?
91  rejectOverlaps = conf.getParameter<bool>("OverlapCleaning");
92 
93  // Split hits ?
94  splitHits = conf.getParameter<bool>("SplitHits");
95 
96  // Reject tracks with several seeds ?
97  // Typically don't do that at HLT for electrons, but do it otherwise
98  seedCleaning = conf.getParameter<bool>("SeedCleaning");
99 
100  simTracks_ = conf.getParameter<edm::InputTag>("SimTracks");
101  estimatorCut_= conf.getParameter<double>("EstimatorCut");
102 }
103 
104 
105 // Virtual destructor needed.
107 
108  if(thePropagator) delete thePropagator;
109 
110  // do nothing
111 #ifdef FAMOS_DEBUG
112  std::cout << "TrackCandidateProducer destructed" << std::endl;
113 #endif
114 
115 
116 }
117 
118 void
120 
121  //services
122 
125 
126  es.get<IdealMagneticFieldRecord>().get(magField);
127  es.get<TrackerDigiGeometryRecord>().get(geometry);
128 
129  theMagField = &(*magField);
130  theGeometry = &(*geometry);
131 
133 }
134 
135  // Functions that gets called by framework every event
136 void
138 
139 #ifdef FAMOS_DEBUG
140  std::cout << "################################################################" << std::endl;
141  std::cout << " TrackCandidateProducer produce init " << std::endl;
142 #endif
143 
144  // Useful typedef's to avoid retyping
145  typedef std::pair<reco::TrackRef,edm::Ref<std::vector<Trajectory> > > TrackPair;
146  typedef std::map<unsigned,TrackPair> TrackMap;
147 
148  // The produced objects
149  std::auto_ptr<TrackCandidateCollection> output(new TrackCandidateCollection);
150  std::auto_ptr<reco::TrackCollection> recoTracks(new reco::TrackCollection);
151  std::auto_ptr<TrackingRecHitCollection> recoHits(new TrackingRecHitCollection);
152  std::auto_ptr<reco::TrackExtraCollection> recoTrackExtras(new reco::TrackExtraCollection);
153  std::auto_ptr<std::vector<Trajectory> > recoTrajectories(new std::vector<Trajectory> );
154  std::auto_ptr<TrajTrackAssociationCollection> recoTrajTrackMap( new TrajTrackAssociationCollection() );
155 
156  // Get the seeds
157  // edm::Handle<TrajectorySeedCollection> theSeeds;
159  e.getByLabel(seedProducer,theSeeds);
160 
161  //Retrieve tracker topology from geometry
163  es.get<IdealGeometryRecord>().get(tTopoHand);
164  const TrackerTopology *tTopo=tTopoHand.product();
165 
166 
167  // No seed -> output an empty track collection
168  if(theSeeds->size() == 0) {
169  e.put(output);
170  e.put(recoTracks);
171  e.put(recoHits);
172  e.put(recoTrackExtras);
173  e.put(recoTrajectories);
174  e.put(recoTrajTrackMap);
175  return;
176  }
177 
178  // Get the GS RecHits
179  // edm::Handle<SiTrackerGSRecHit2DCollection> theGSRecHits;
181  e.getByLabel(hitProducer, theGSRecHits);
182 
183  //get other general things
184  const std::vector<unsigned> theSimTrackIds = theGSRecHits->ids();
185  // SimTracks and SimVertices
186 
188  e.getByLabel("famosSimHits",theSimVtx);
190  e.getByLabel("famosSimHits",theSTC);
191 
192  const edm::SimTrackContainer* theSimTracks = &(*theSTC);
193  LogDebug("FastTracking")<<"looking at: "<< theSimTrackIds.size()<<" simtracks.";
194 
195 
196 
197 
198  // The input track collection + extra's
199  /*
200  edm::Handle<reco::TrackCollection> theTrackCollection;
201  edm:: Handle<std::vector<Trajectory> > theTrajectoryCollection;
202  edm::Handle<TrajTrackAssociationCollection> theAssoMap;
203  bool isTrackCollection = e.getByLabel(trackProducer,theTrackCollection);
204  */
205  std::vector<edm::Handle<reco::TrackCollection> > theTrackCollections;
206  std::vector<edm:: Handle<std::vector<Trajectory> > > theTrajectoryCollections;
207  std::vector<edm::Handle<TrajTrackAssociationCollection> > theAssoMaps;
208  std::vector<bool> isTrackCollections;
211  TrackMap theTrackMap;
212  unsigned nCollections = trackProducers.size();
213  unsigned nRecoHits = 0;
214 
215  if ( nCollections ) {
216  theTrackCollections.resize(nCollections);
217  theTrajectoryCollections.resize(nCollections);
218  theAssoMaps.resize(nCollections);
219  isTrackCollections.resize(nCollections);
220  for ( unsigned tprod=0; tprod < nCollections; ++tprod ) {
221  isTrackCollections[tprod] = e.getByLabel(trackProducers[tprod],theTrackCollections[tprod]);
222 
223  if ( isTrackCollections[tprod] ) {
224  // The track collection
225  reco::TrackCollection::const_iterator aTrack = theTrackCollections[tprod]->begin();
226  reco::TrackCollection::const_iterator lastTrack = theTrackCollections[tprod]->end();
227  // The numbers of hits
228  for ( ; aTrack!=lastTrack; ++aTrack ) nRecoHits+= aTrack->recHitsSize();
229  e.getByLabel(trackProducers[tprod],theTrajectoryCollections[tprod]);
230  e.getByLabel(trackProducers[tprod],theAssoMaps[tprod]);
231  // The association between trajectories and tracks
232  anAssociation = theAssoMaps[tprod]->begin();
233  lastAssociation = theAssoMaps[tprod]->end();
234 #ifdef FAMOS_DEBUG
235  std::cout << "Input Track Producer " << tprod << " : " << trackProducers[tprod] << std::endl;
236  std::cout << "List of tracks already reconstructed " << std::endl;
237 #endif
238  // Build the map of correspondance between reco tracks and sim tracks
239  for ( ; anAssociation != lastAssociation; ++anAssociation ) {
240  edm::Ref<std::vector<Trajectory> > aTrajectoryRef = anAssociation->key;
241  reco::TrackRef aTrackRef = anAssociation->val;
242  // Find the simtrack id of the reconstructed track
243  int recoTrackId = findId(*aTrackRef);
244  if ( recoTrackId < 0 ) continue;
245 #ifdef FAMOS_DEBUG
246  std::cout << recoTrackId << " ";
247 #endif
248  // And store it.
249  theTrackMap[recoTrackId] = TrackPair(aTrackRef,aTrajectoryRef);
250  }
251 #ifdef FAMOS_DEBUG
252  std::cout << std::endl;
253 #endif
254  }
255  }
256  // This is to save some time at push_back.
257  recoHits->reserve(nRecoHits);
258  }
259 
260 
261  // Loop over the seeds
262  int currentTrackId = -1;
263  /*
264  TrajectorySeedCollection::const_iterator aSeed = theSeeds->begin();
265  TrajectorySeedCollection::const_iterator lastSeed = theSeeds->end();
266  for ( ; aSeed!=lastSeed; ++aSeed ) {
267  // The first hit of the seed and its simtrack id
268  */
269  /* */
270 #ifdef FAMOS_DEBUG
271  std::cout << "Input seed Producer : " << seedProducer << std::endl;
272  std::cout << "Number of seeds : " << theSeeds->size() << std::endl;
273 #endif
274  unsigned seed_size = theSeeds->size();
275  for (unsigned seednr = 0; seednr < seed_size; ++seednr){
276 
277  LogDebug("FastTracking")<<"looking at seed #:"<<seednr;
278 
279  // The seed
280  const BasicTrajectorySeed* aSeed = &((*theSeeds)[seednr]);
281 
282  std::vector<int> simTrackIds;
283  std::map<int,TrajectoryStateOnSurface> seedStates;
284  std::map<int,TrajectoryStateOnSurface> simtkStates;
285 
286  TrackerRecHit theFirstSeedingTrackerRecHit;
287  if (theSeeds->at(seednr).nHits()==0){
288  //new stuff for no hits on seed
289 
290  LogDebug("FastTracking")<<" seed with no hits to be considered.";
291 
292  //moved out of the loop
293  //edm::ESHandle<MagneticField> field;
294  //es.get<IdealMagneticFieldRecord>().get(field);
295 
296  PTrajectoryStateOnDet ptod =theSeeds->at(seednr).startingState();
297  DetId id(ptod.detId());
298  const GeomDet * g = theGeometry->idToDet(id);
299  const Surface * surface=&g->surface();
300 
302 
304  es.get<TrackingComponentsRecord>().get("AnyDirectionAnalyticalPropagator",propagator);
305 
306  //moved out of the loop
307  // const std::vector<unsigned> theSimTrackIds = theGSRecHits->ids();
308  // edm::Handle<edm::SimTrackContainer> theSTC;
309  // e.getByLabel(simTracks_,theSTC);
310  // const edm::SimTrackContainer* theSimTracks = &(*theSTC);
311 
312  double minimunEst=1000000;
313  LogDebug("FastTracking")<<"looking at: "<< theSimTrackIds.size()<<" simtracks.";
314  for ( unsigned tkId=0; tkId != theSimTrackIds.size(); ++tkId ) {
315 
316  const SimTrack & simtrack = (*theSimTracks)[theSimTrackIds[tkId]];
317 
319  simtrack.trackerSurfacePosition().y(),
320  simtrack.trackerSurfacePosition().z());
321 
322  GlobalVector momentum(simtrack.trackerSurfaceMomentum().x(),
323  simtrack.trackerSurfaceMomentum().y(),
324  simtrack.trackerSurfaceMomentum().z());
325 
326  if (position.basicVector().dot( momentum.basicVector() ) * seedState.globalPosition().basicVector().dot(seedState.globalMomentum().basicVector()) <0. ){
327  LogDebug("FastTracking")<<"not on the same direction.";
328  continue;
329  }
330 
331  //no charge mis-identification ... FIXME
332  int charge = (int) simtrack.charge();
333  GlobalTrajectoryParameters glb_parameters(position,
334  momentum,
335  charge,
336  theMagField);
337  FreeTrajectoryState simtrack_trackerstate(glb_parameters);
338 
339  TrajectoryStateOnSurface simtrack_comparestate = propagator->propagate(simtrack_trackerstate,*surface);
340 
341 
342  if (!simtrack_comparestate.isValid()){
343  LogDebug("FastTracking")<<" ok this is a state-based seed. simtrack state does not propagate to the seed surface. skipping.";
344  continue;}
345 
346  if (simtrack_comparestate.globalPosition().basicVector().dot(simtrack_comparestate.globalMomentum().basicVector()) * seedState.globalPosition().basicVector().dot(seedState.globalMomentum().basicVector()) <0. ){
347  LogDebug("FastTracking")<<"not on the same direction.";
348  continue;
349  }
350 
351 
352  AlgebraicVector5 v(seedState.localParameters().vector() - simtrack_comparestate.localParameters().vector());
353  AlgebraicSymMatrix55 m(seedState.localError().matrix());
354  bool ierr = !m.Invert();
355  if ( ierr ){
356  edm::LogWarning("FastTracking") <<" Candidate Producer cannot invert the error matrix! - Skipping...";
357  continue;
358  }
359  double est = ROOT::Math::Similarity(v,m);
360  LogDebug("FastTracking")<<"comparing two state on the seed surface:\n"
361  <<"seed: "<<seedState
362  <<"sim: "<<simtrack_comparestate
363  <<"\n estimator is: "<<est;
364 
365  if (est<minimunEst) minimunEst=est;
366  if (est<estimatorCut_){
367  simTrackIds.push_back(theSimTrackIds[tkId]);
368  //making a state with exactly the sim track parameters
369  //the initial errors are set to unity just for kicks
370  // AlgebraicSymMatrix C(5,1);// C*=50;
371  //new attempt!!!!
373  C *= 0.0000001;
374 
375  seedStates[theSimTrackIds[tkId]] = TrajectoryStateOnSurface(simtrack_comparestate.globalParameters(),
377  seedState.surface());
378  LogDebug("FastTracking")<<"the compatibility estimator is: "<<est<<" for track id: "<<simTrackIds.back();
379  }
380  }//SimTrack loop
381  if (simTrackIds.size()==0) LogDebug("FastTracking")<<"could not find any simtrack within errors, closest was at: "<<minimunEst;
382  }//seed has 0 hit.
383  else{
384  //same old stuff
385  // Find the first hit of the Seed
386  TrajectorySeed::range theSeedingRecHitRange = aSeed->recHits();
387  const SiTrackerGSMatchedRecHit2D * theFirstSeedingRecHit = (const SiTrackerGSMatchedRecHit2D*) (&(*(theSeedingRecHitRange.first)));
388  theFirstSeedingTrackerRecHit = TrackerRecHit(theFirstSeedingRecHit,theGeometry,tTopo);
389  // The SimTrack id associated to that recHit
390  simTrackIds.push_back( theFirstSeedingRecHit->simtrackId() );
391  }
392 
393  //from then on, only the simtrack IDs are usefull.
394  //now loop over all possible trackid for this seed.
395  //an actual seed can be shared by two tracks in dense envirronement, and also for hit-less seeds.
396  for (unsigned int iToMake=0;iToMake!=simTrackIds.size();++iToMake){
397  int simTrackId = simTrackIds[iToMake];
398 
399  // Don't consider seeds belonging to a track already considered
400  // (Equivalent to seed cleaning)
401  if ( seedCleaning && currentTrackId == simTrackId ) continue;
402  currentTrackId = simTrackId;
403 
404  // A vector of TrackerRecHits belonging to the track and the number of crossed layers
405  std::vector<TrackerRecHit> theTrackerRecHits;
406  unsigned theNumberOfCrossedLayers = 0;
407 
408  // The track has indeed been reconstructed already -> Save the pertaining info
409  TrackMap::const_iterator theTrackIt = theTrackMap.find(simTrackId);
410  if ( nCollections && theTrackIt != theTrackMap.end() ) {
411 
412  if ( keepFittedTracks ) {
413  LogDebug("FastTracking") << "Track " << simTrackId << " already reconstructed -> copy it";
414  // The track and trajectroy references
415  reco::TrackRef aTrackRef = theTrackIt->second.first;
416  edm::Ref<std::vector<Trajectory> > aTrajectoryRef = theTrackIt->second.second;
417 
418  // A copy of the track
419  reco::Track aRecoTrack(*aTrackRef);
420  recoTracks->push_back(aRecoTrack);
421 
422  // A copy of the hits
423  unsigned nh = aRecoTrack.recHitsSize();
424  for ( unsigned ih=0; ih<nh; ++ih ) {
425  TrackingRecHit *hit = aRecoTrack.recHit(ih)->clone();
426  recoHits->push_back(hit);
427  }
428 
429  // A copy of the trajectories
430  recoTrajectories->push_back(*aTrajectoryRef);
431 
432  }// keepFitterTracks
433  else {
434  LogDebug("FastTracking") << "Track " << simTrackId << " already reconstructed -> ignore it";
435  }
436 
437  // The track was not saved -> create a track candidate.
438 
439  } //already existing collection of tracks
440  else{//no collection of tracks already exists
441 
442  LogDebug("FastTracking")<<"Track " << simTrackId << " is considered to return a track candidate" ;
443 
444  // Get all the rechits associated to this track
445  SiTrackerGSMatchedRecHit2DCollection::range theRecHitRange = theGSRecHits->get(simTrackId);
446  SiTrackerGSMatchedRecHit2DCollection::const_iterator theRecHitRangeIteratorBegin = theRecHitRange.first;
447  SiTrackerGSMatchedRecHit2DCollection::const_iterator theRecHitRangeIteratorEnd = theRecHitRange.second;
449 
450  LogDebug("FastTracking")<<"counting: "<<theRecHitRangeIteratorEnd-theRecHitRangeIteratorBegin<<" hits to be considered.";
451 
452  bool firstRecHit = true;
453  TrackerRecHit theCurrentRecHit, thePreviousRecHit;
454  TrackerRecHit theFirstHitComp, theSecondHitComp;
455 
456  for ( iterRecHit = theRecHitRangeIteratorBegin;
457  iterRecHit != theRecHitRangeIteratorEnd;
458  ++iterRecHit) {
459 
460  // Check the number of crossed layers
461  if ( theNumberOfCrossedLayers >= maxNumberOfCrossedLayers ) continue;
462 
463  // Get current and previous rechits
464  if(!firstRecHit) thePreviousRecHit = theCurrentRecHit;
465  theCurrentRecHit = TrackerRecHit(&(*iterRecHit),theGeometry,tTopo);
466 
467  //>>>>>>>>>BACKBUILDING CHANGE: DO NOT STAT FROM THE FIRST HIT OF THE SEED
468 
469  // NOTE: checking the direction --> specific for OIHit only
470  // if( aSeed->direction()!=oppositeToMomentum ) {
471  // // Check that the first rechit is indeed the first seeding hit
472  // if ( firstRecHit && theCurrentRecHit != theFirstSeedingTrackerRecHit && theSeeds->at(seednr).nHits()!=0 ) continue;
473  // }
474  //>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
475 
476  // Count the number of crossed layers
477  if ( !theCurrentRecHit.isOnTheSameLayer(thePreviousRecHit) )
478  ++theNumberOfCrossedLayers;
479 
480  // Add all rechits (Grouped Trajectory Builder) from this hit onwards
481  // Always add the first seeding rechit anyway
482  if ( !rejectOverlaps || firstRecHit ) {
483  // Split matched hits (if requested / possible )
484  if ( splitHits && theCurrentRecHit.matchedHit()->isMatched() ) addSplitHits(theCurrentRecHit,theTrackerRecHits);
485  else theTrackerRecHits.push_back(theCurrentRecHit); // No splitting
486  firstRecHit = false;
487 
488  // And now treat the following RecHits if hits in the same layer
489  // have to be rejected - The split option is not
490  } else {
491 
492  // Not the same layer : Add the current hit
493  if ( theCurrentRecHit.subDetId() != thePreviousRecHit.subDetId() ||
494  theCurrentRecHit.layerNumber() != thePreviousRecHit.layerNumber() ) {
495 
496  // Split matched hits (if requested / possible )
497  if ( splitHits && theCurrentRecHit.matchedHit()->isMatched() ) addSplitHits(theCurrentRecHit,theTrackerRecHits);
498  else theTrackerRecHits.push_back(theCurrentRecHit); // No splitting
499 
500  // Same layer : keep the current hit if better, and drop the other - otherwise do nothing
501  } else if ( theCurrentRecHit.localError() < thePreviousRecHit.localError() ) {
502 
503  // Split matched hits (if requested / possible )
504  if( splitHits && theCurrentRecHit.matchedHit()->isMatched() ){
505  // Remove the previous hit(s)
506  theTrackerRecHits.pop_back();
507  if ( thePreviousRecHit.matchedHit()->isMatched() ) theTrackerRecHits.pop_back();
508 
509  // Replace by the new hits
510  addSplitHits(theCurrentRecHit,theTrackerRecHits);
511  }
512  // No splitting
513  else {
514  theTrackerRecHits.back() = theCurrentRecHit; // Replace the previous hit by the current hit
515  }
516 
517  } else {
518 
519  //keep the older rechit as a reference if the error of the new one is worse
520  theCurrentRecHit = thePreviousRecHit;
521  }
522  }
523  }
524 
525  // End of loop over the track rechits
526  }//no collection of track already existed. adding the hits by hand.
527 
528  LogDebug("FastTracking")<<" number of hits: " << theTrackerRecHits.size()<<" after counting overlaps and splitting.";
529 
530  // 1) Create the OwnWector of TrackingRecHits
532  unsigned nTrackerHits = theTrackerRecHits.size();
533  recHits.reserve(nTrackerHits); // To save some time at push_back
534 
535  if (aSeed->direction()==oppositeToMomentum){
536  LogDebug("FastTracking")<<"reversing the order of the hits";
537  std::reverse(theTrackerRecHits.begin(),theTrackerRecHits.end());
538  }
539 
540  for ( unsigned ih=0; ih<nTrackerHits; ++ih ) {
541  TrackingRecHit* aTrackingRecHit = theTrackerRecHits[ih].hit()->clone();
542  recHits.push_back(aTrackingRecHit);
543 
544  const DetId& detId = theTrackerRecHits[ih].hit()->geographicalId();
545  LogDebug("FastTracking")
546  << "Added RecHit from detid " << detId.rawId()
547  << " subdet = " << theTrackerRecHits[ih].subDetId()
548  << " layer = " << theTrackerRecHits[ih].layerNumber()
549  << " ring = " << theTrackerRecHits[ih].ringNumber()
550  << " error = " << theTrackerRecHits[ih].localError()
551  << std::endl
552 
553  << "Track/z/r : "
554  << simTrackId << " "
555  << theTrackerRecHits[ih].globalPosition().z() << " "
556  << theTrackerRecHits[ih].globalPosition().perp() << std::endl;
557  if ( theTrackerRecHits[ih].matchedHit() && theTrackerRecHits[ih].matchedHit()->isMatched() )
558  LogTrace("FastTracking") << "Matched : " << theTrackerRecHits[ih].matchedHit()->isMatched()
559  << "Rphi Hit = " << theTrackerRecHits[ih].matchedHit()->monoHit()->simhitId()
560  << "Stereo Hit = " << theTrackerRecHits[ih].matchedHit()->stereoHit()->simhitId()
561  <<std::endl;
562  }//loop over the rechits
563 
564  // Check the number of crossed layers
565  if ( theNumberOfCrossedLayers < minNumberOfCrossedLayers ) {
566  LogDebug("FastTracking")<<"not enough layer crossed ("<<theNumberOfCrossedLayers<<")";
567  continue;
568  }
569 
570 
571  //>>>>>>>>>BACKBUILDING CHANGE: REPLACE THE STARTING STATE
572 
573  // Create a track Candidate (now with the reference to the seed!) .
574  //PTrajectoryStateOnDet PTSOD = aSeed->startingState();
575  PTrajectoryStateOnDet PTSOD;
576 
577  if (aSeed->nHits()==0){
578  //stabilize the fit with the true simtrack state
579  //in case of zero hits
580 
581  PTSOD = trajectoryStateTransform::persistentState(seedStates[simTrackId],aSeed->startingState().detId());
582 
583  } else {
584  //create the initial state from the SimTrack
585  int vertexIndex = (*theSimTracks)[currentTrackId].vertIndex();
586  // a) origin vertex
587  GlobalPoint position((*theSimVtx)[vertexIndex].position().x(),
588  (*theSimVtx)[vertexIndex].position().y(),
589  (*theSimVtx)[vertexIndex].position().z());
590 
591  // b) initial momentum
592  GlobalVector momentum( (*theSimTracks)[currentTrackId].momentum().x() ,
593  (*theSimTracks)[currentTrackId].momentum().y() ,
594  (*theSimTracks)[currentTrackId].momentum().z() );
595  // c) electric charge
596  float charge = (*theSimTracks)[simTrackId].charge();
597  // -> inital parameters
598  GlobalTrajectoryParameters initialParams(position,momentum,(int)charge,theMagField);
599  // -> large initial errors
600  AlgebraicSymMatrix55 errorMatrix= AlgebraicMatrixID();
601  CurvilinearTrajectoryError initialError(errorMatrix);
602  // -> initial state
603  FreeTrajectoryState initialFTS(initialParams, initialError);
604 #ifdef FAMOS_DEBUG
605  std::cout << "TrajectorySeedProducer: FTS momentum " << initialFTS.momentum() << std::endl;
606 #endif
607  const GeomDet* initialLayer = theGeometry->idToDet(recHits.front().geographicalId());
608  //this is wrong because the FTS is defined at vertex, and it need to be properly propagated to the first rechit
609  // const TrajectoryStateOnSurface initialTSOS(initialFTS, initialLayer->surface());
610  const TrajectoryStateOnSurface initialTSOS = thePropagator->propagate(initialFTS,initialLayer->surface()) ;
611  if (!initialTSOS.isValid()) continue;
612 
613 
614  PTSOD = trajectoryStateTransform::persistentState(initialTSOS,recHits.front().geographicalId().rawId());
615  }
616 
618  newTrackCandidate(recHits,
619  *aSeed,
620  PTSOD,
621  edm::RefToBase<TrajectorySeed>(theSeeds,seednr));
622 
623  LogDebug("FastTracking")<< "\tSeed Information " << std::endl
624  << "\tSeed Direction = " << aSeed->direction() << std::endl
625  << "\tSeed StartingDet = " << aSeed->startingState().detId() << std::endl
626  << "\tTrajectory Parameters " << std::endl
627  << "\t\t detId = " << newTrackCandidate.trajectoryStateOnDet().detId() << std::endl
628  << "\t\t loc.px = "
629  << newTrackCandidate.trajectoryStateOnDet().parameters().momentum().x()
630  << std::endl
631  << "\t\t loc.py = "
632  << newTrackCandidate.trajectoryStateOnDet().parameters().momentum().y()
633  << std::endl
634  << "\t\t loc.pz = "
635  << newTrackCandidate.trajectoryStateOnDet().parameters().momentum().z()
636  << std::endl
637  << "\t\t error = ";
638 
639  bool newTrackCandidateIsDuplicate = isDuplicateCandidate(*output,newTrackCandidate);
640  if (!newTrackCandidateIsDuplicate) output->push_back(newTrackCandidate);
641  LogDebug("FastTracking")<<"filling a track candidate into the collection, now having: "<<output->size();
642 
643  }//loop over possible simtrack associated.
644  }//loop over all possible seeds.
645 
646  // Save the track candidates in the event
647  LogDebug("FastTracking") << "Saving "
648  << output->size() << " track candidates and "
649  << recoTracks->size() << " reco::Tracks ";
650  // Save the track candidates
651  e.put(output);
652 
653 
654 
655  // Save the tracking recHits
656 
657  edm::OrphanHandle <TrackingRecHitCollection> theRecoHits = e.put(recoHits );
658 
659  // Create the track extras and add the references to the rechits
660  unsigned hits=0;
661  unsigned nTracks = recoTracks->size();
662  recoTrackExtras->reserve(nTracks); // To save some time at push_back
663  for ( unsigned index = 0; index < nTracks; ++index ) {
664  //reco::TrackExtra aTrackExtra;
665  reco::Track& aTrack = recoTracks->at(index);
666  reco::TrackExtra aTrackExtra(aTrack.outerPosition(),
667  aTrack.outerMomentum(),
668  aTrack.outerOk(),
669  aTrack.innerPosition(),
670  aTrack.innerMomentum(),
671  aTrack.innerOk(),
672  aTrack.outerStateCovariance(),
673  aTrack.outerDetId(),
674  aTrack.innerStateCovariance(),
675  aTrack.innerDetId(),
676  aTrack.seedDirection(),
677  aTrack.seedRef());
678 
679  unsigned nHits = aTrack.recHitsSize();
680  for ( unsigned int ih=0; ih<nHits; ++ih) {
681  aTrackExtra.add(TrackingRecHitRef(theRecoHits,hits++));
682  }
683  recoTrackExtras->push_back(aTrackExtra);
684  }
685 
686 
687  // Save the track extras
688  edm::OrphanHandle<reco::TrackExtraCollection> theRecoTrackExtras = e.put(recoTrackExtras);
689 
690  // Add the reference to the track extra in the tracks
691  for ( unsigned index = 0; index<nTracks; ++index ) {
692  const reco::TrackExtraRef theTrackExtraRef(theRecoTrackExtras,index);
693  (recoTracks->at(index)).setExtra(theTrackExtraRef);
694  }
695 
696  // Save the tracks
697  edm::OrphanHandle<reco::TrackCollection> theRecoTracks = e.put(recoTracks);
698 
699  // Save the trajectories
700  edm::OrphanHandle<std::vector<Trajectory> > theRecoTrajectories = e.put(recoTrajectories);
701 
702  // Create and set the trajectory/track association map
703  for ( unsigned index = 0; index<nTracks; ++index ) {
704  edm::Ref<std::vector<Trajectory> > trajRef( theRecoTrajectories, index );
705  edm::Ref<reco::TrackCollection> tkRef( theRecoTracks, index );
706  recoTrajTrackMap->insert(trajRef,tkRef);
707  }
708 
709 
710  // Save the association map.
711  e.put(recoTrajTrackMap);
712 
713 }
714 
715 int
717  int trackId = -1;
718  trackingRecHit_iterator aHit = aTrack.recHitsBegin();
719  trackingRecHit_iterator lastHit = aTrack.recHitsEnd();
720  for ( ; aHit!=lastHit; ++aHit ) {
721  if ( !aHit->get()->isValid() ) continue;
722  // const SiTrackerGSRecHit2D * rechit = (const SiTrackerGSRecHit2D*) (aHit->get());
723  const SiTrackerGSMatchedRecHit2D * rechit = (const SiTrackerGSMatchedRecHit2D*) (aHit->get());
724  trackId = rechit->simtrackId();
725  break;
726  }
727  return trackId;
728 }
729 
730 void
732  std::vector<TrackerRecHit>& theTrackerRecHits) {
733 
734  const SiTrackerGSRecHit2D* mHit = theCurrentRecHit.matchedHit()->monoHit();
735  const SiTrackerGSRecHit2D* sHit = theCurrentRecHit.matchedHit()->stereoHit();
736 
737  // Add the new hits
738  if( mHit->simhitId() < sHit->simhitId() ) {
739 
740  theTrackerRecHits.push_back(TrackerRecHit(mHit,theCurrentRecHit));
741  theTrackerRecHits.push_back(TrackerRecHit(sHit,theCurrentRecHit));
742 
743  } else {
744 
745  theTrackerRecHits.push_back(TrackerRecHit(sHit,theCurrentRecHit));
746  theTrackerRecHits.push_back(TrackerRecHit(mHit,theCurrentRecHit));
747 
748  }
749 
750 }
751 
753  typedef TrackCandidateCollection::const_iterator TCCI;
754 
755  TCCI candsEnd = candidates.end();
756  double newQbp = newCand.trajectoryStateOnDet().parameters().qbp();
757  double newDxdz = newCand.trajectoryStateOnDet().parameters().dxdz();
758  TrackCandidate::range newHits = newCand.recHits();
759 
760  for (TCCI iCand = candidates.begin(); iCand!= candsEnd; ++iCand){
761  //pick some first traits of duplication: qbp and dxdz
762  double iQbp = iCand->trajectoryStateOnDet().parameters().qbp();
763  double iDxdz = iCand->trajectoryStateOnDet().parameters().dxdz();
764  if (newQbp == iQbp && newDxdz == iDxdz){
765  LogDebug("isDuplicateCandidate")<<"Probably a duplicate "<<iQbp <<" : "<<iDxdz;
766  TrackCandidate::range iHits = iCand->recHits();
767  unsigned int nHits = 0;
768  unsigned int nShared = 0;
769  unsigned int nnewHits = 0;
770  for (TrackCandidate::const_iterator iiHit = iHits.first; iiHit != iHits.second; ++iiHit){
771  nHits++;
772  for (TrackCandidate::const_iterator inewHit = newHits.first; inewHit != newHits.second; ++inewHit){
773  if (iiHit == iHits.first) nnewHits++;
774  if (sameLocalParameters(&*iiHit,&*inewHit)) nShared++;
775  }
776  }
777  LogDebug("isDuplicateCandidate") <<"nHits "<<nHits<<" nnewHits "<< nnewHits<< " nShared "<<nShared;
778  if (nHits == nShared && nHits == nnewHits) return true;
779  }
780  }
781  return false;
782 }
783 
785  bool localPointSame = aH->localPosition() == bH->localPosition();
786  bool localErrorSame = aH->localPositionError().xx() == bH->localPositionError().xx();
787  localErrorSame &= aH->localPositionError().yy() == bH->localPositionError().yy();
788  localErrorSame &= aH->localPositionError().xy() == bH->localPositionError().xy();
789  return localPointSame && localErrorSame;
790 }
#define LogDebug(id)
PropagationDirection direction() const
const math::XYZVectorD & trackerSurfacePosition() const
Definition: SimTrack.h:36
T getParameter(std::string const &) const
float xx() const
Definition: LocalError.h:24
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:50
const int & simhitId() const
const LocalTrajectoryParameters & localParameters() const
range recHits() const
size_t recHitsSize() const
Get number of RecHits. (Warning, this includes invalid hits, which are not physical hits)...
Definition: Track.h:68
ROOT::Math::SMatrixIdentity AlgebraicMatrixID
TrackCandidateProducer(const edm::ParameterSet &conf)
int findId(const reco::Track &aTrack) const
const CurvilinearTrajectoryError & curvilinearError() const
std::vector< TrackCandidate > TrackCandidateCollection
size_type size() const
Definition: OwnVector.h:247
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:10
PTrajectoryStateOnDet persistentState(const TrajectoryStateOnSurface &ts, unsigned int detid)
T y() const
Definition: PV3DBase.h:63
GlobalPoint globalPosition() const
bool innerOk() const
return true if the innermost hit is valid
Definition: Track.h:39
const SiTrackerGSRecHit2D * stereoHit() const
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
std::pair< const_iterator, const_iterator > range
const math::XYZPoint & outerPosition() const
position of the outermost hit
Definition: Track.h:46
float charge() const
charge
Definition: CoreSimTrack.cc:3
const MagneticField * theMagField
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:35
double charge(const std::vector< uint8_t > &Ampls)
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
float float float z
static int position[TOTALCHAMBERS][3]
Definition: ReadPGInfo.cc:509
bool isOnTheSameLayer(const TrackerRecHit &other) const
Check if two hits are on the same layer of the same subdetector.
virtual void produce(edm::Event &e, const edm::EventSetup &es) override
const math::XYZPoint & innerPosition() const
position of the innermost hit
Definition: Track.h:41
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
PTrajectoryStateOnDet const & trajectoryStateOnDet() const
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
AlgebraicVector5 vector() const
void push_back(D *&d)
Definition: OwnVector.h:273
float xy() const
Definition: LocalError.h:25
const SurfaceType & surface() const
const SiTrackerGSMatchedRecHit2D * matchedHit() const
The Hit itself.
Definition: TrackerRecHit.h:73
const SiTrackerGSRecHit2D * monoHit() const
float yy() const
Definition: LocalError.h:26
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
PropagatorWithMaterial * thePropagator
edm::Ref< TrackingRecHitCollection > TrackingRecHitRef
persistent reference to a TrackingRecHit
CovarianceMatrix outerStateCovariance() const
outermost trajectory state curvilinear errors
Definition: Track.h:50
T z() const
Definition: PV3DBase.h:64
std::vector< edm::InputTag > trackProducers
unsigned int outerDetId() const
DetId of the detector on which surface the outermost state is located.
Definition: Track.h:58
std::pair< const_iterator, const_iterator > range
bool isMatched(TrackingRecHit const &hit)
unsigned int detId() const
const AlgebraicSymMatrix55 & matrix() const
LocalVector momentum() const
Momentum vector in the local frame.
const LocalTrajectoryError & localError() const
trackingRecHit_iterator recHitsBegin() const
Iterator to first hit on the track.
Definition: Track.h:62
edm::AssociationMap< edm::OneToOne< std::vector< Trajectory >, reco::TrackCollection, unsigned short > > TrajTrackAssociationCollection
virtual const GeomDet * idToDet(DetId) const
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:390
GlobalVector momentum() const
#define LogTrace(id)
tuple conf
Definition: dbtoconf.py:185
virtual TrackingRecHit * clone() const =0
Definition: DetId.h:18
std::vector< TrackExtra > TrackExtraCollection
collection of TrackExtra objects
Definition: TrackExtraFwd.h:9
edm::RefToBase< TrajectorySeed > seedRef() const
Definition: Track.h:111
PTrajectoryStateOnDet const & startingState() const
const math::XYZVector & outerMomentum() const
momentum vector at the outermost hit position
Definition: Track.h:48
TrajectoryStateOnSurface transientState(const PTrajectoryStateOnDet &ts, const Surface *surface, const MagneticField *field)
ROOT::Math::SVector< double, 5 > AlgebraicVector5
bool outerOk() const
return true if the outermost hit is valid
Definition: Track.h:37
const GlobalTrajectoryParameters & globalParameters() const
virtual LocalError localPositionError() const =0
const math::XYZTLorentzVectorD & trackerSurfaceMomentum() const
Definition: SimTrack.h:38
virtual void beginRun(edm::Run const &run, const edm::EventSetup &es) override
const T & get() const
Definition: EventSetup.h:55
CovarianceMatrix innerStateCovariance() const
innermost trajectory state curvilinear errors
Definition: Track.h:52
key_type key() const
Accessor for product key.
Definition: Ref.h:266
T const * product() const
Definition: ESHandle.h:62
range recHits() const
unsigned int layerNumber() const
The Layer Number.
Definition: TrackerRecHit.h:84
ESHandle< TrackerGeometry > geometry
unsigned int nHits() const
const AlgebraicSymMatrix55 & matrix() const
GlobalVector globalMomentum() const
const math::XYZVector & innerMomentum() const
momentum vector at the innermost hit position
Definition: Track.h:44
TrackingRecHitRef recHit(size_t i) const
Get i-th hit on the track.
Definition: Track.h:66
unsigned int subDetId() const
The subdet Id.
Definition: TrackerRecHit.h:81
double localError()
PropagationDirection seedDirection() const
direction of how the hits were sorted in the original seed
Definition: Track.h:104
tuple cout
Definition: gather_cfg.py:121
const TrackerGeometry * theGeometry
DetId geographicalId() const
Definition: DDAxes.h:10
bool sameLocalParameters(const TrackingRecHit *aH, const TrackingRecHit *bH) const
bool isDuplicateCandidate(const TrackCandidateCollection &candidates, const TrackCandidate &newCand) const
RecHitContainer::const_iterator const_iterator
T x() const
Definition: PV3DBase.h:62
virtual LocalPoint localPosition() const =0
const BasicVectorType & basicVector() const
Definition: PV3DBase.h:56
unsigned int innerDetId() const
DetId of the detector on which surface the innermost state is located.
Definition: Track.h:60
std::vector< SimTrack > SimTrackContainer
reference front()
Definition: OwnVector.h:348
void reserve(size_t)
Definition: OwnVector.h:267
const LocalTrajectoryParameters & parameters() const
Definition: Run.h:41
void addSplitHits(const TrackerRecHit &, std::vector< TrackerRecHit > &)
T dot(const Basic3DVector &rh) const
Scalar product, or &quot;dot&quot; product, with a vector of same type.
trackingRecHit_iterator recHitsEnd() const
Iterator to last hit on the track.
Definition: Track.h:64