CMS 3D CMS Logo

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 
00015 /*****************************************************************************/
00016 TrackListCombiner::TrackListCombiner(const edm::ParameterSet& ps)
00017 {
00018   trackProducers = ps.getParameter<vector<string> >("trackProducers");
00019 
00020   produces<reco::TrackCollection>();
00021   produces<vector<Trajectory> >();
00022 }
00023 
00024 /*****************************************************************************/
00025 TrackListCombiner::~TrackListCombiner()
00026 {
00027 }
00028 
00029 /*****************************************************************************/
00030 void TrackListCombiner::produce(edm::Event& ev, const edm::EventSetup& es)
00031 {
00032   auto_ptr<reco::TrackCollection> outputTracks(new reco::TrackCollection);
00033   auto_ptr<vector<Trajectory> >   outputTrajes(new vector<Trajectory>   );
00034 
00035   LogTrace("MinBiasTracking")
00036     << "[TrackListCombiner]";
00037 
00038   // Go through all track producers
00039   int i = 1;
00040   for(vector<string>::iterator trackProducer = trackProducers.begin();
00041                                trackProducer!= trackProducers.end();
00042                                trackProducer++, i++)
00043   {
00044     reco::TrackBase::TrackAlgorithm algo;
00045     switch(i) 
00046     {
00047       case 1:  algo = reco::TrackBase::iter1; break;
00048       case 2:  algo = reco::TrackBase::iter2; break;
00049       case 3:  algo = reco::TrackBase::iter3; break;
00050       default: algo = reco::TrackBase::undefAlgorithm;
00051     }
00052 
00053     // Get track collection
00054     edm::Handle<reco::TrackCollection> trackHandle;
00055     ev.getByLabel(*trackProducer,      trackHandle);
00056     const reco::TrackCollection & trackCollection = *(trackHandle.product());
00057 
00058     // Get trajectory collection
00059     edm::Handle<vector<Trajectory> > trajeHandle;
00060     ev.getByLabel(*trackProducer,    trajeHandle);
00061     const vector<Trajectory> & trajeCollection = *(trajeHandle.product());
00062 
00063     LogTrace("MinBiasTracking")
00064       << " [TrackListCombiner] " << *trackProducer
00065       << " : " << trackCollection.size() << "|" << trajeCollection.size();
00066 
00067     for(reco::TrackCollection::const_iterator track = trackCollection.begin();
00068                                               track!= trackCollection.end();
00069                                               track++)
00070     {
00071       // Get new track
00072       reco::Track * theTrack = new reco::Track(*track);
00073 
00074       // Set extra
00075       reco::TrackExtraRef theTrackExtraRef = track->extra();    
00076       theTrack->setExtra(theTrackExtraRef);    
00077       theTrack->setHitPattern((*theTrackExtraRef).recHits());
00078      
00079       // Set algorithm
00080       theTrack->setAlgorithm(algo);
00081 
00082       // Store track
00083       outputTracks->push_back(*theTrack);
00084 
00085       delete theTrack;
00086     }
00087 
00088     for(vector<Trajectory>::const_iterator traje = trajeCollection.begin();
00089                                            traje!= trajeCollection.end();
00090                                            traje++)
00091     {
00092       // Get new trajectory
00093       Trajectory * theTraje = new Trajectory(*traje);
00094 
00095       // Store trajectory
00096       outputTrajes->push_back(*theTraje);
00097 
00098       delete theTraje;
00099     }
00100   }
00101 
00102   LogTrace("MinBiasTracking")
00103     << " [TrackListCombiner] allTracks : " << outputTracks->size()
00104                                     << "|" << outputTrajes->size();
00105 
00106   // Put back result to event
00107   ev.put(outputTracks);
00108   ev.put(outputTrajes);
00109 }

Generated on Tue Jun 9 17:44:51 2009 for CMSSW by  doxygen 1.5.4