CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/RecoPixelVertexing/PixelLowPtUtilities/plugins/TrackListCombiner.cc

Go to the documentation of this file.
00001 #include "TrackListCombiner.h"
00002 
00003 #include "FWCore/Framework/interface/Event.h"
00004 #include "FWCore/Framework/interface/EventSetup.h"
00005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00006 
00007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00008 
00009 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00010 #include "DataFormats/TrackReco/interface/Track.h"
00011 #include "DataFormats/TrackReco/interface/TrackExtra.h"
00012 
00013 #include "TrackingTools/PatternTools/interface/Trajectory.h"
00014 #include "TrackingTools/PatternTools/interface/TrajTrackAssociation.h"
00015 
00016 using namespace std;
00017 
00018 /*****************************************************************************/
00019 TrackListCombiner::TrackListCombiner(const edm::ParameterSet& ps)
00020 {
00021   trackProducers = ps.getParameter<vector<string> >("trackProducers");
00022 
00023   produces<reco::TrackCollection>();
00024   produces<reco::TrackExtraCollection>();
00025   produces<TrackingRecHitCollection>();
00026   produces<vector<Trajectory> >();
00027   produces<TrajTrackAssociationCollection>();
00028 }
00029 
00030 /*****************************************************************************/
00031 TrackListCombiner::~TrackListCombiner()
00032 {
00033 }
00034 
00035 /*****************************************************************************/
00036 void TrackListCombiner::produce(edm::Event& ev, const edm::EventSetup& es)
00037 {
00038   auto_ptr<reco::TrackCollection>          recoTracks
00039       (new reco::TrackCollection);
00040   auto_ptr<reco::TrackExtraCollection>     recoTrackExtras
00041       (new reco::TrackExtraCollection);
00042   auto_ptr<TrackingRecHitCollection>       recoHits
00043       (new TrackingRecHitCollection);
00044   auto_ptr<vector<Trajectory> >            recoTrajectories
00045       (new vector<Trajectory>);
00046   auto_ptr<TrajTrackAssociationCollection> recoTrajTrackMap
00047       (new TrajTrackAssociationCollection());
00048 
00049   LogTrace("MinBiasTracking")
00050     << "[TrackListCombiner]";
00051 
00052   // Go through all track producers
00053   int i = 1;
00054   for(vector<string>::iterator trackProducer = trackProducers.begin();
00055                                trackProducer!= trackProducers.end();
00056                                trackProducer++, i++)
00057   {
00058     reco::TrackBase::TrackAlgorithm algo;
00059     switch(i) 
00060     {
00061       case 1:  algo = reco::TrackBase::iter1; break;
00062       case 2:  algo = reco::TrackBase::iter2; break;
00063       case 3:  algo = reco::TrackBase::iter3; break;
00064       default: algo = reco::TrackBase::undefAlgorithm;
00065     }
00066 
00067     edm::Handle<vector<Trajectory> > theTrajectoryCollection;
00068     edm::Handle<TrajTrackAssociationCollection> theAssoMap;  
00069 
00070     ev.getByLabel(*trackProducer, theTrajectoryCollection);
00071     ev.getByLabel(*trackProducer, theAssoMap);
00072 
00073     LogTrace("MinBiasTracking")
00074       << " [TrackListCombiner] " << *trackProducer
00075       << " : " << theAssoMap->size();
00076 
00077     
00078     // The track collection iterators
00079     TrajTrackAssociationCollection::const_iterator anAssociation;  
00080     TrajTrackAssociationCollection::const_iterator lastAssociation;
00081     anAssociation = theAssoMap->begin();
00082     lastAssociation = theAssoMap->end();
00083 
00084     // Build the map of correspondance between reco tracks and sim tracks
00085     for ( ; anAssociation != lastAssociation; ++anAssociation )
00086     { 
00087       edm::Ref<vector<Trajectory> > aTrajectoryRef = anAssociation->key;
00088       reco::TrackRef aTrackRef = anAssociation->val;
00089       
00090       // A copy of the track
00091       reco::Track aRecoTrack(*aTrackRef);
00092 
00093       // Set algorithm
00094       aRecoTrack.setAlgorithm(algo);
00095 
00096       recoTracks->push_back(aRecoTrack);      
00097 
00098       // A copy of the hits
00099       unsigned nh = aRecoTrack.recHitsSize();
00100       for(unsigned ih=0; ih<nh; ++ih)
00101       {
00102         TrackingRecHit *hit = aRecoTrack.recHit(ih)->clone();
00103         recoHits->push_back(hit);
00104       }
00105       
00106       // A copy of the trajectories
00107       recoTrajectories->push_back(*aTrajectoryRef);
00108       
00109     }
00110   }
00111 
00112   LogTrace("MinBiasTracking")
00113     << " [TrackListCombiner] allTracks : " << recoTracks->size()
00114                                     << "|" << recoTrajectories->size();
00115 
00116   // Save the tracking recHits
00117   edm::OrphanHandle<TrackingRecHitCollection> theRecoHits = ev.put(recoHits);
00118   
00119   // Create the track extras and add the references to the rechits
00120   unsigned hits = 0;
00121   unsigned nTracks = recoTracks->size();
00122   recoTrackExtras->reserve(nTracks); // To save some time at push_back
00123   for(unsigned index = 0; index < nTracks; ++index )
00124   { 
00125     reco::Track& aTrack = recoTracks->at(index);
00126     reco::TrackExtra aTrackExtra(aTrack.outerPosition(),
00127                                  aTrack.outerMomentum(),
00128                                  aTrack.outerOk(),
00129                                  aTrack.innerPosition(),
00130                                  aTrack.innerMomentum(),
00131                                  aTrack.innerOk(),
00132                                  aTrack.outerStateCovariance(),
00133                                  aTrack.outerDetId(),
00134                                  aTrack.innerStateCovariance(),
00135                                  aTrack.innerDetId(),
00136                                  aTrack.seedDirection(),
00137                                  aTrack.seedRef());
00138     
00139     unsigned nHits = aTrack.recHitsSize();
00140     for ( unsigned int ih=0; ih<nHits; ++ih)
00141       aTrackExtra.add(TrackingRecHitRef(theRecoHits,hits++));
00142     recoTrackExtras->push_back(aTrackExtra);
00143   }
00144   
00145   // Save the track extras
00146   edm::OrphanHandle<reco::TrackExtraCollection> theRecoTrackExtras =
00147     ev.put(recoTrackExtras);
00148   
00149   // Add the reference to the track extra in the tracks
00150   for(unsigned index = 0; index<nTracks; ++index)
00151   { 
00152     const reco::TrackExtraRef theTrackExtraRef(theRecoTrackExtras,index);
00153     (recoTracks->at(index)).setExtra(theTrackExtraRef);
00154   }
00155   
00156   // Save the tracks
00157   edm::OrphanHandle<reco::TrackCollection> theRecoTracks = ev.put(recoTracks);
00158   
00159   // Save the trajectories
00160   edm::OrphanHandle<vector<Trajectory> > theRecoTrajectories =
00161     ev.put(recoTrajectories);
00162   
00163   // Create and set the trajectory/track association map 
00164   for(unsigned index = 0; index<nTracks; ++index)
00165   { 
00166     edm::Ref<vector<Trajectory> > trajRef( theRecoTrajectories, index );
00167     edm::Ref<reco::TrackCollection>    tkRef( theRecoTracks, index );
00168     recoTrajTrackMap->insert(trajRef,tkRef);
00169   }
00170   
00171   // Save the association map
00172   ev.put(recoTrajTrackMap);
00173 }
00174