CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch12/src/FastSimulation/Tracking/plugins/FastTrackMerger.cc

Go to the documentation of this file.
00001 #include <memory>
00002 
00003 #include "FWCore/Framework/interface/Event.h"
00004 #include "FWCore/Framework/interface/EventSetup.h"
00005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00006 #include "FWCore/Framework/interface/ESHandle.h"
00007 
00008 #include "DataFormats/Common/interface/Handle.h"
00009 #include "DataFormats/Common/interface/OwnVector.h"
00010 #include "DataFormats/TrackerRecHit2D/interface/SiTrackerGSRecHit2DCollection.h" 
00011 #include "DataFormats/TrackerRecHit2D/interface/SiTrackerGSMatchedRecHit2DCollection.h" 
00012 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00013 #include "DataFormats/TrackReco/interface/Track.h"
00014 #include "DataFormats/TrackReco/interface/TrackExtra.h"
00015 #include "DataFormats/TrackReco/interface/TrackExtraFwd.h"
00016 
00017 #include "TrackingTools/PatternTools/interface/Trajectory.h"
00018 #include "TrackingTools/PatternTools/interface/TrajTrackAssociation.h"
00019 
00020 #include "FastSimulation/Tracking/plugins/FastTrackMerger.h"
00021 
00022 #include <vector>
00023 #include <map>
00024 //
00025 
00026 //for debug only 
00027 //#define FAMOS_DEBUG
00028 
00029 FastTrackMerger::FastTrackMerger(const edm::ParameterSet& conf) 
00030 {  
00031 #ifdef FAMOS_DEBUG
00032   std::cout << "FastTrackMerger created" << std::endl;
00033 #endif
00034 
00035   // The main product is a track collection, and all extras
00036   produces<reco::TrackCollection>();
00037   
00038   // The name of the track producers to merge
00039   trackProducers = conf.getParameter<std::vector<edm::InputTag> >("TrackProducers");
00040 
00041   // The name of the track producers to remove
00042   std::vector<edm::InputTag> defaultRemove;
00043   removeTrackProducers = 
00044     conf.getUntrackedParameter<std::vector<edm::InputTag> >("RemoveTrackProducers",defaultRemove);
00045 
00046   // Only the tracks!
00047   tracksOnly = conf.getUntrackedParameter<bool>("SaveTracksOnly",false);
00048 
00049   // optional pT cut
00050   double pTMin = conf.getUntrackedParameter<bool>("pTMin",0.);
00051   pTMin2 = pTMin*pTMin;
00052 
00053   // optional nHit cut
00054   minHits = conf.getUntrackedParameter<unsigned>("minHits",0);
00055 
00056   // optional track quality saving
00057   promoteQuality = conf.getUntrackedParameter<bool>("promoteTrackQuality",false);
00058   qualityStr = conf.getUntrackedParameter<std::string>("newQuality","");
00059 
00060   // optional trackAlgo (iter0/1/2/3/4)
00061   trackAlgo = conf.getUntrackedParameter<unsigned>("trackAlgo",0);
00062 
00063   //new parameters for Trajectory filtering
00064 
00065   // The minimum number of hits
00066   theMinimumNumberOfHits = conf.getUntrackedParameter<unsigned>("MinNumberOfTrajHits",0);
00067 
00068   // The maximum number of Lost Hits
00069   theMaxLostHits = conf.getUntrackedParameter<unsigned>("MaxLostTrajHits",99);
00070 
00071   // the max number of consecutive Lost Hits 
00072   theMaxConsecutiveLostHits = conf.getUntrackedParameter<unsigned>("MaxConsecutiveLostTrajHits",3);
00073 
00074   //======================
00075 
00076   if ( !tracksOnly ) { 
00077     produces<reco::TrackExtraCollection>();
00078     produces<TrackingRecHitCollection>();
00079     produces<std::vector<Trajectory> >();
00080     produces<TrajTrackAssociationCollection>();
00081   } else {
00082     produces<reco::TrackExtraCollection>();
00083     produces<TrackingRecHitCollection>();
00084   }
00085 }
00086 
00087   
00088 // Functions that gets called by framework every event
00089 void 
00090 FastTrackMerger::produce(edm::Event& e, const edm::EventSetup& es) {        
00091 
00092 #ifdef FAMOS_DEBUG
00093   std::cout << "################################################################" << std::endl;
00094   std::cout << " FastTrackMerger produce init " << std::endl;
00095 #endif
00096 
00097   // The track algorithm (only of merging after iterative tracking)
00098   reco::TrackBase::TrackAlgorithm algo;
00099   switch(trackAlgo)  {
00100   case 4:  algo = reco::TrackBase::iter0; break;
00101   case 5:  algo = reco::TrackBase::iter1; break;
00102   case 6:  algo = reco::TrackBase::iter2; break;
00103   case 7:  algo = reco::TrackBase::iter3; break;
00104   case 8:  algo = reco::TrackBase::iter4; break;
00105   case 9:  algo = reco::TrackBase::iter5; break;
00106   default: algo = reco::TrackBase::undefAlgorithm;
00107   }
00108 
00109   // The produced objects
00110   std::auto_ptr<reco::TrackCollection> recoTracks(new reco::TrackCollection);
00111   std::auto_ptr<reco::TrackExtraCollection> recoTrackExtras(new reco::TrackExtraCollection);
00112   std::auto_ptr<TrackingRecHitCollection> recoHits(new TrackingRecHitCollection);
00113   std::auto_ptr<std::vector<Trajectory> > recoTrajectories(new std::vector<Trajectory>);
00114   std::auto_ptr<TrajTrackAssociationCollection> recoTrajTrackMap(new TrajTrackAssociationCollection());
00115 
00116   //prepare the ref to copy the extras
00117   reco::TrackExtraRefProd rTrackExtras = e.getRefBeforePut<reco::TrackExtraCollection>();
00118   TrackingRecHitRefProd rHits = e.getRefBeforePut<TrackingRecHitCollection>();
00119   std::vector<reco::TrackRef> trackRefs;
00120 
00121   // No seed -> output an empty track collection
00122   if(trackProducers.size() == 0) {
00123     e.put(recoTracks);
00124     if ( !tracksOnly ) { 
00125       e.put(recoTrackExtras);
00126       e.put(recoHits);
00127       e.put(recoTrajectories);
00128       e.put(recoTrajTrackMap);
00129     } else {
00130       e.put(recoTrackExtras);
00131       e.put(recoHits);
00132     }
00133     return;
00134   }
00135 
00136   // The quality to be set
00137   reco::TrackBase::TrackQuality qualityToSet;
00138   if (qualityStr != "")
00139     qualityToSet = reco::TrackBase::qualityByName(qualityStr);
00140   else 
00141     qualityToSet = reco::TrackBase::undefQuality;
00142 
00143   // The input track collection handle
00144   edm::Handle<reco::TrackCollection> theTrackCollection;
00145 
00146   // First, the tracks to be removed
00147   std::set<unsigned> removeTracks;
00148   for ( unsigned aProducer=0; aProducer<removeTrackProducers.size(); ++aProducer ) { 
00149     bool isTrackCollection = e.getByLabel(removeTrackProducers[aProducer],theTrackCollection); 
00150     if (!isTrackCollection) continue;
00151     reco::TrackCollection::const_iterator aTrack = theTrackCollection->begin();
00152     reco::TrackCollection::const_iterator lastTrack = theTrackCollection->end();
00153     for ( ; aTrack!=lastTrack; ++aTrack ) {
00154       // Get the simtrack Id
00155       int recoTrackId = findId(*aTrack);
00156       if ( recoTrackId < 0 ) continue;
00157       // Remove the track if requested
00158       if ( removeTracks.find((unsigned)recoTrackId) != removeTracks.end() ) continue;
00159       removeTracks.insert((unsigned)recoTrackId);
00160     }      
00161   }
00162   
00163   // Then the tracks to be added
00164   std::set<unsigned> alreadyAddedTracks;
00165   
00166   // Loop on the track producers to be merged
00167   for ( unsigned aProducer=0; aProducer<trackProducers.size(); ++aProducer ) { 
00168     
00169     bool isTrackCollection = e.getByLabel(trackProducers[aProducer],theTrackCollection); 
00170     if ( ! isTrackCollection ) { 
00171 #ifdef FAMOS_DEBUG
00172       std::cout << "***FastTrackMerger*** Warning! The track collection " 
00173                 << trackProducers[aProducer].encode() 
00174                 << " does not exist." << std::endl;
00175 #endif
00176       continue;
00177     }
00178 
00179 #ifdef FAMOS_DEBUG
00180     std::cout << "***FastTrackMerger*** of track collection " 
00181               << trackProducers[aProducer].encode() 
00182               << " with " << theTrackCollection->size() 
00183               << " tracks to be copied"
00184               << std::endl;
00185 #endif
00186     reco::TrackCollection::const_iterator aTrack = theTrackCollection->begin();
00187     reco::TrackCollection::const_iterator lastTrack = theTrackCollection->end();
00188 
00189     //Copy Tracks and TracksExtra 
00190     if ( tracksOnly ) { 
00191 
00192       edm:: Handle<reco::TrackExtraCollection > theTrackExtraCollection;
00193       //      bool isTrackExtraCollection = e.getByLabel(trackProducers[aProducer],theTrackExtraCollection); 
00194 
00195       bool index = 0;
00196       for ( ; aTrack!=lastTrack; ++aTrack,++index ) {
00197 
00198         // Find the track id
00199         int recoTrackId = findId(*aTrack);
00200         if ( recoTrackId < 0 ) continue;
00201         
00202         // Ignore tracks to be removed or tracks already copied
00203         std::set<unsigned>::iterator iR = removeTracks.find((unsigned)recoTrackId);
00204 #ifdef FAMOS_DEBUG
00205         if( iR != removeTracks.end() ) std::cout << recoTrackId << "(REMOVED) ";
00206 #endif
00207         if( iR != removeTracks.end() ) continue;
00208         
00209         // Ignore tracks already copied
00210         std::set<unsigned>::iterator iA = alreadyAddedTracks.find((unsigned)recoTrackId);
00211 #ifdef FAMOS_DEBUG
00212         if( iA != alreadyAddedTracks.end() ) std::cout << recoTrackId << "(ALREADY ADDED) ";
00213 #endif
00214         if( iA != alreadyAddedTracks.end() ) continue;
00215         //if it is not there then add it! 
00216         alreadyAddedTracks.insert(recoTrackId);
00217 
00218         
00219 #ifdef FAMOS_DEBUG
00220         std::cout << recoTrackId << " ADDED ";
00221 #endif
00222 
00223         // Ignore tracks with too small a pT
00224 #ifdef FAMOS_DEBUG
00225         if ( aTrack->innerMomentum().Perp2() < pTMin2 ) 
00226           std::cout << "PTMIN CUT APPLIED = " <<  aTrack->innerMomentum().Perp2() << std::endl;
00227 #endif
00228 
00229         if ( aTrack->innerMomentum().Perp2() < pTMin2 ) continue;
00230         
00231         // Ignore tracks with too small a pT
00232         if ( aTrack->recHitsSize() < minHits ) continue;
00233         
00234         // A copy of the track + save the transient reference to the track extra reference
00235         reco::Track aRecoTrack(*aTrack);
00236         recoTracks->push_back(aRecoTrack);
00237         
00238 
00239         recoTrackExtras->push_back(reco::TrackExtra(aRecoTrack.outerPosition(), aRecoTrack.outerMomentum(), aRecoTrack.outerOk(),
00240                                                aRecoTrack.innerPosition(), aRecoTrack.innerMomentum(), aRecoTrack.innerOk(),
00241                                                aRecoTrack.outerStateCovariance(), aRecoTrack.outerDetId(),
00242                                                aRecoTrack.innerStateCovariance(), aRecoTrack.innerDetId(),
00243                                                aRecoTrack.seedDirection(), aRecoTrack.seedRef() ) ); 
00244 
00245 
00246         recoTracks->back().setExtra(reco::TrackExtraRef(rTrackExtras,recoTrackExtras->size()-1)); 
00247 
00248         reco::TrackExtra & tx = recoTrackExtras->back();
00249         tx.setResiduals(aRecoTrack.residuals());
00250         //TrackingRecHits
00251         for( trackingRecHit_iterator hit = aRecoTrack.recHitsBegin(); hit != aRecoTrack.recHitsEnd(); ++ hit ) {
00252           recoHits->push_back( (*hit)->clone() );
00253           tx.add( TrackingRecHitRef( rHits , recoHits->size() - 1) );
00254         }
00255         
00256         // Save the quality if requested
00257         if (promoteQuality) recoTracks->back().setQuality(qualityToSet);
00258         if ( trackAlgo ) recoTracks->back().setAlgorithm(algo);
00259         
00260       }
00261       
00262       
00263       // All extras are to be copied too -> loop on the Trajectory/Track map association 
00264     } else { 
00265       
00266       edm:: Handle<std::vector<Trajectory> > theTrajectoryCollection;
00267       edm::Handle<TrajTrackAssociationCollection> theAssoMap;  
00268 
00269       // Count the number of hits and reserve appropriate memory
00270       unsigned nRecoHits = 0;
00271       for ( ; aTrack!=lastTrack; ++aTrack ) nRecoHits+= aTrack->recHitsSize();
00272       recoHits->reserve(nRecoHits); // This is to save some time at push_back.
00273       
00274       e.getByLabel(trackProducers[aProducer],theTrajectoryCollection);
00275       e.getByLabel(trackProducers[aProducer],theAssoMap);
00276       
00277       // The track collection iterators.
00278       TrajTrackAssociationCollection::const_iterator anAssociation;  
00279       TrajTrackAssociationCollection::const_iterator lastAssociation;
00280       anAssociation = theAssoMap->begin();
00281       lastAssociation = theAssoMap->end();
00282 #ifdef FAMOS_DEBUG
00283       std::cout << "List of tracks to be copied " << std::endl;
00284 #endif
00285       // Build the map of correspondance between reco tracks and sim tracks
00286       for ( ; anAssociation != lastAssociation; ++anAssociation ) { 
00287         edm::Ref<std::vector<Trajectory> > aTrajectoryRef = anAssociation->key;
00288         reco::TrackRef aTrackRef = anAssociation->val;
00289         // Find the track id
00290         int recoTrackId = findId(*aTrackRef);
00291         if ( recoTrackId < 0 ) continue;
00292         
00293         // Ignore tracks to be removed or tracks already copied
00294         std::set<unsigned>::iterator iR = removeTracks.find((unsigned)recoTrackId);
00295 #ifdef FAMOS_DEBUG
00296         if( iR != removeTracks.end() ) std::cout << recoTrackId << "(REMOVED) ";
00297 #endif
00298         if( iR != removeTracks.end() ) continue;
00299         
00300         // Ignore tracks already copied
00301         std::set<unsigned>::iterator iA = alreadyAddedTracks.find((unsigned)recoTrackId);
00302 #ifdef FAMOS_DEBUG
00303         if( iA != alreadyAddedTracks.end() ) std::cout << recoTrackId << "(ALREADY ADDED) ";
00304 #endif
00305         if( iA != alreadyAddedTracks.end() ) continue;
00306         
00307 #ifdef FAMOS_DEBUG
00308         std::cout << recoTrackId << " Newly Added " << std::endl;
00309 #endif
00310         //if it is not there then add it! 
00311         alreadyAddedTracks.insert(recoTrackId);
00312         
00313         // Ignore tracks with too small a pT
00314         if ( aTrackRef->innerMomentum().Perp2() < pTMin2 ) continue;
00315         
00316         // Ignore tracks with too few hits
00317         if ( aTrackRef->recHitsSize() < minHits ) continue;
00318         
00319         //==== add more cuts on the trajectory (emulate the Trajectory Filter) 
00320 
00321 #ifdef FAMOS_DEBUG
00322         if(aTrajectoryRef->lostHits() > theMaxLostHits )
00323           std::cout << "\tmaxLostHits= " << aTrajectoryRef->lostHits() << "\tCUT =" << theMaxLostHits << std::endl;
00324 #endif
00325 
00326         if((unsigned)aTrajectoryRef->lostHits() > theMaxLostHits ) continue;
00327 
00328 #ifdef FAMOS_DEBUG
00329         if(aTrajectoryRef->foundHits() < theMinimumNumberOfHits )
00330           std::cout << "\tMinimumNumberOfHits = " <<  aTrajectoryRef->foundHits() << "\tCUT = " <<theMinimumNumberOfHits <<  std::endl;
00331 #endif
00332         
00333         if((unsigned)aTrajectoryRef->foundHits() < theMinimumNumberOfHits ) continue;
00334         //calculate the consecutive Lost Hits
00335         unsigned consecLostHits = 0;
00336         const std::vector<TrajectoryMeasurement> tms = aTrajectoryRef->measurements();
00337         for(int itm= tms.size();itm!=0; --itm){
00338           if(tms[itm-1].recHit()->isValid())break;
00339           else if (Trajectory::lost(*tms[itm-1].recHit())) consecLostHits++;
00340         }
00341 
00342 #ifdef FAMOS_DEBUG
00343         if( consecLostHits > theMaxConsecutiveLostHits ) 
00344           std::cout << "\tconsecLostHits = " << consecLostHits << std::endl;
00345 #endif
00346 
00347         if( consecLostHits > theMaxConsecutiveLostHits ) continue;
00348 
00349 
00350         //=============end new filters
00351 
00352 
00353         // A copy of the track
00354         reco::Track aRecoTrack(*aTrackRef);
00355         recoTracks->push_back(aRecoTrack);      
00356         // Save the quality if requested
00357         if (promoteQuality) recoTracks->back().setQuality(qualityToSet);        
00358         if ( trackAlgo ) recoTracks->back().setAlgorithm(algo);
00359         
00360      
00361         // A copy of the hits
00362         unsigned nh = aRecoTrack.recHitsSize();
00363         for ( unsigned ih=0; ih<nh; ++ih ) {
00364           TrackingRecHit *hit = aRecoTrack.recHit(ih)->clone();
00365           recoHits->push_back(hit);
00366         }
00367         
00368         // A copy of the trajectories
00369         recoTrajectories->push_back(*aTrajectoryRef);
00370         
00371       }
00372 #ifdef FAMOS_DEBUG
00373       std::cout << std::endl;
00374 #endif
00375     }
00376 }
00377     
00378     // Save the track candidates in the event
00379 #ifdef FAMOS_DEBUG
00380   std::cout << "Saving " 
00381             << recoTracks->size() << " reco::Tracks " << std::endl;
00382 
00383 #endif
00384   
00385   if ( tracksOnly ) { 
00386     // Save only the tracks (with transient reference to track extras)
00387     e.put(recoTracks);
00388     e.put(recoTrackExtras);
00389     e.put(recoHits);
00390 
00391 
00392   } else { 
00393     // Save the tracking recHits
00394     edm::OrphanHandle <TrackingRecHitCollection> theRecoHits = e.put(recoHits );
00395     
00396     // Create the track extras and add the references to the rechits
00397     unsigned hits=0;
00398     unsigned nTracks = recoTracks->size();
00399     recoTrackExtras->reserve(nTracks); // To save some time at push_back
00400     for ( unsigned index = 0; index < nTracks; ++index ) { 
00401 
00402       //reco::TrackExtra aTrackExtra;
00403       reco::Track& aTrack = recoTracks->at(index);
00404 
00405       // std::cout << "initial TrackAlgo = " << trackAlgo << "\tNAME " << aTrack.algoName() << "\tnumber = " << aTrack.algo() << std::endl;
00406 
00407       reco::TrackExtra aTrackExtra(aTrack.outerPosition(),
00408                                    aTrack.outerMomentum(),
00409                                    aTrack.outerOk(),
00410                                    aTrack.innerPosition(),
00411                                    aTrack.innerMomentum(),
00412                                    aTrack.innerOk(),
00413                                    aTrack.outerStateCovariance(),
00414                                    aTrack.outerDetId(),
00415                                    aTrack.innerStateCovariance(),
00416                                    aTrack.innerDetId(),
00417                                    aTrack.seedDirection(),
00418                                    aTrack.seedRef());
00419       
00420       unsigned nHits = aTrack.recHitsSize();
00421       for ( unsigned int ih=0; ih<nHits; ++ih) {
00422         aTrackExtra.add(TrackingRecHitRef(theRecoHits,hits++));
00423       }
00424       recoTrackExtras->push_back(aTrackExtra);
00425     }
00426     
00427     // Save the track extras
00428     edm::OrphanHandle<reco::TrackExtraCollection> theRecoTrackExtras = e.put(recoTrackExtras);
00429     
00430     // Add the reference to the track extra in the tracks
00431     for ( unsigned index = 0; index<nTracks; ++index ) { 
00432       const reco::TrackExtraRef theTrackExtraRef(theRecoTrackExtras,index);
00433       (recoTracks->at(index)).setExtra(theTrackExtraRef);
00434     }
00435     
00436     // Save the tracks
00437     edm::OrphanHandle<reco::TrackCollection> theRecoTracks = e.put(recoTracks);
00438     
00439     // Save the trajectories
00440     edm::OrphanHandle<std::vector<Trajectory> > theRecoTrajectories = e.put(recoTrajectories);
00441     
00442     // Create and set the trajectory/track association map 
00443     for ( unsigned index = 0; index<nTracks; ++index ) { 
00444       edm::Ref<std::vector<Trajectory> > trajRef( theRecoTrajectories, index );
00445       edm::Ref<reco::TrackCollection>    tkRef( theRecoTracks, index );
00446       recoTrajTrackMap->insert(trajRef,tkRef);
00447     }
00448     
00449     // Save the association map.
00450     e.put(recoTrajTrackMap);
00451     
00452   }
00453 
00454 }
00455 
00456 int 
00457 FastTrackMerger::findId(const reco::Track& aTrack) const {
00458   int trackId = -1;
00459   trackingRecHit_iterator aHit = aTrack.recHitsBegin();
00460   trackingRecHit_iterator lastHit = aTrack.recHitsEnd();
00461   for ( ; aHit!=lastHit; ++aHit ) {
00462     if ( !aHit->get()->isValid() ) continue;
00463     //    const SiTrackerGSRecHit2D * rechit = (const SiTrackerGSRecHit2D*) (aHit->get());
00464     const SiTrackerGSMatchedRecHit2D * rechit = (const SiTrackerGSMatchedRecHit2D*) (aHit->get());
00465     trackId = rechit->simtrackId();
00466     break;
00467   }
00468   return trackId;
00469 }
00470 
00471