CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/FastSimulation/Tracking/plugins/PixelTracksProducer.cc

Go to the documentation of this file.
00001 #include <memory>
00002 #include "FastSimulation/Tracking/plugins/PixelTracksProducer.h"
00003 
00004 #include "FWCore/Framework/interface/Event.h"
00005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00006 
00007 #include "DataFormats/Common/interface/Handle.h"
00008 #include "DataFormats/Common/interface/OwnVector.h"
00009 
00010 #include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"
00011 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
00012 
00013 //Pixel Specific stuff
00014 #include "RecoTracker/TkTrackingRegions/interface/TrackingRegionProducer.h"
00015 #include "RecoTracker/TkTrackingRegions/interface/TrackingRegionProducerFactory.h"
00016 
00017 #include "RecoPixelVertexing/PixelTrackFitting/interface/PixelFitter.h"
00018 #include "RecoPixelVertexing/PixelTrackFitting/interface/PixelFitterFactory.h"
00019 
00020 #include "RecoPixelVertexing/PixelTrackFitting/interface/PixelTrackFilter.h"
00021 #include "RecoPixelVertexing/PixelTrackFitting/interface/PixelTrackFilterFactory.h"
00022 #include "RecoPixelVertexing/PixelTrackFitting/interface/TracksWithHits.h"
00023 
00024 #include "DataFormats/TrackReco/interface/Track.h"
00025 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00026 #include "DataFormats/TrackReco/interface/TrackExtra.h"
00027 #include "DataFormats/Common/interface/OrphanHandle.h"
00028 #include <vector>
00029 
00030 using namespace pixeltrackfitting;
00031 
00032 PixelTracksProducer::PixelTracksProducer(const edm::ParameterSet& conf) : 
00033   theFitter(0), 
00034   theFilter(0), 
00035   theRegionProducer(0)
00036 {  
00037 
00038   produces<reco::TrackCollection>();
00039   produces<TrackingRecHitCollection>();
00040   produces<reco::TrackExtraCollection>();
00041 
00042   const edm::ParameterSet& regfactoryPSet = conf.getParameter<edm::ParameterSet>("RegionFactoryPSet");
00043   std::string regfactoryName = regfactoryPSet.getParameter<std::string>("ComponentName");
00044   theRegionProducer = TrackingRegionProducerFactory::get()->create(regfactoryName,regfactoryPSet);
00045   
00046   const edm::ParameterSet& fitterPSet = conf.getParameter<edm::ParameterSet>("FitterPSet");
00047   std::string fitterName = fitterPSet.getParameter<std::string>("ComponentName");
00048   theFitter = PixelFitterFactory::get()->create( fitterName, fitterPSet);
00049   
00050   const edm::ParameterSet& filterPSet = conf.getParameter<edm::ParameterSet>("FilterPSet");
00051   std::string filterName = filterPSet.getParameter<std::string>("ComponentName");
00052   theFilter = PixelTrackFilterFactory::get()->create( filterName, filterPSet);
00053   
00054   // The name of the seed producer
00055   seedProducer = conf.getParameter<edm::InputTag>("SeedProducer");
00056 
00057 }
00058 
00059   
00060 // Virtual destructor needed.
00061 PixelTracksProducer::~PixelTracksProducer() {
00062 
00063   delete theFilter;
00064   delete theFitter;
00065   delete theRegionProducer;
00066 
00067 } 
00068  
00069 
00070 // Functions that gets called by framework every event
00071 void 
00072 PixelTracksProducer::produce(edm::Event& e, const edm::EventSetup& es) {        
00073   
00074   std::auto_ptr<reco::TrackCollection> tracks(new reco::TrackCollection);    
00075   std::auto_ptr<TrackingRecHitCollection> recHits(new TrackingRecHitCollection);
00076   std::auto_ptr<reco::TrackExtraCollection> trackExtras(new reco::TrackExtraCollection);
00077   typedef std::vector<const TrackingRecHit *> RecHits;
00078   
00079   TracksWithRecHits pixeltracks;
00080   TracksWithRecHits cleanedTracks;
00081   
00082   edm::Handle<TrajectorySeedCollection> theSeeds;
00083   e.getByLabel(seedProducer,theSeeds);
00084 
00085   // No seed -> output an empty track collection
00086   if(theSeeds->size() == 0) {
00087     e.put(tracks);
00088     e.put(recHits);
00089     e.put(trackExtras);
00090     return;
00091   }
00092   
00093   //only one region Global, but it is called at every event...
00094   //maybe there is a smarter way to set it only once
00095   //NEED TO FIX
00096   typedef std::vector<TrackingRegion* > Regions;
00097   typedef Regions::const_iterator IR;
00098   Regions regions = theRegionProducer->regions(e,es);
00099   for (IR ir=regions.begin(), irEnd=regions.end(); ir < irEnd; ++ir) {
00100     const TrackingRegion & region = **ir;
00101     
00102     // Loop over the seeds
00103     TrajectorySeedCollection::const_iterator aSeed = theSeeds->begin();
00104     TrajectorySeedCollection::const_iterator lastSeed = theSeeds->end();
00105     for ( ; aSeed!=lastSeed; ++aSeed ) { 
00106       
00107       // Find the first hit and last hit of the Seed
00108       TrajectorySeed::range theSeedingRecHitRange = aSeed->recHits();
00109       edm::OwnVector<TrackingRecHit>::const_iterator aSeedingRecHit = theSeedingRecHitRange.first;
00110       edm::OwnVector<TrackingRecHit>::const_iterator theLastSeedingRecHit = theSeedingRecHitRange.second;
00111 
00112       // Loop over the rechits
00113       std::vector<const TrackingRecHit*> TripletHits(3,static_cast<const TrackingRecHit*>(0));
00114       for ( unsigned i=0; aSeedingRecHit!=theLastSeedingRecHit; ++i,++aSeedingRecHit )  
00115         TripletHits[i] = &(*aSeedingRecHit);
00116       
00117       // fitting the triplet
00118       reco::Track* track = theFitter->run(es, TripletHits, region);
00119       
00120       // decide if track should be skipped according to filter 
00121       if ( ! (*theFilter)(track) ) { 
00122         delete track; 
00123         continue; 
00124       }
00125       
00126       // add tracks 
00127       pixeltracks.push_back(TrackWithRecHits(track, TripletHits));
00128       
00129     }
00130   }
00131   
00132   int cc=0;
00133   int nTracks = pixeltracks.size();
00134   for (int i = 0; i < nTracks; ++i) {
00135 
00136     reco::Track* track   =  pixeltracks.at(i).first;
00137     const RecHits & hits = pixeltracks.at(i).second;
00138     
00139     for (unsigned int k = 0; k < hits.size(); k++) {
00140       TrackingRecHit *hit = (hits.at(k))->clone();
00141       track->setHitPattern(*hit, k);
00142       recHits->push_back(hit);
00143     }
00144 
00145     tracks->push_back(*track);
00146     delete track;
00147     
00148   }
00149   
00150   edm::OrphanHandle <TrackingRecHitCollection> ohRH = e.put( recHits );
00151   
00152   for (int k = 0; k < nTracks; ++k) {
00153 
00154     // reco::TrackExtra* theTrackExtra = new reco::TrackExtra();
00155     reco::TrackExtra theTrackExtra;
00156     
00157     //fill the TrackExtra with TrackingRecHitRef
00158     // unsigned int nHits = tracks->at(k).numberOfValidHits();
00159     unsigned nHits = 3; // We are dealing with triplets!
00160     for(unsigned int i = 0; i < nHits; ++i) {
00161       theTrackExtra.add(TrackingRecHitRef(ohRH,cc++));
00162       //theTrackExtra->add(TrackingRecHitRef(ohRH,cc));
00163       //cc++;
00164     }
00165     
00166     trackExtras->push_back(theTrackExtra);
00167     //trackExtras->push_back(*theTrackExtra);
00168     //delete theTrackExtra;
00169   }
00170   
00171   edm::OrphanHandle<reco::TrackExtraCollection> ohTE = e.put(trackExtras);
00172   
00173   for (int k = 0; k < nTracks; k++) {
00174 
00175     const reco::TrackExtraRef theTrackExtraRef(ohTE,k);
00176     (tracks->at(k)).setExtra(theTrackExtraRef);
00177 
00178   }
00179   
00180   e.put(tracks);
00181   
00182   // Avoid a memory leak !
00183   unsigned nRegions = regions.size();
00184   for ( unsigned iRegions=0; iRegions<nRegions; ++iRegions ) {
00185     delete regions[iRegions];
00186   }
00187 
00188 }
00189