CMS 3D CMS Logo

PixelTrackProducer.cc

Go to the documentation of this file.
00001 #include "PixelTrackProducer.h"
00002 
00003 #include "FWCore/Framework/interface/Event.h"
00004 #include "FWCore/Framework/interface/EventSetup.h"
00005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00006 
00007 
00008 #include "RecoTracker/TkTrackingRegions/interface/OrderedHitsGenerator.h"
00009 #include "RecoTracker/TkTrackingRegions/interface/OrderedHitsGeneratorFactory.h"
00010 
00011 #include "RecoTracker/TkTrackingRegions/interface/TrackingRegion.h"
00012 #include "RecoTracker/TkTrackingRegions/interface/TrackingRegionProducer.h"
00013 #include "RecoTracker/TkTrackingRegions/interface/TrackingRegionProducerFactory.h"
00014 
00015 #include "RecoPixelVertexing/PixelTrackFitting/interface/PixelFitter.h"
00016 #include "RecoPixelVertexing/PixelTrackFitting/interface/PixelFitterFactory.h"
00017 
00018 #include "RecoPixelVertexing/PixelTrackFitting/interface/PixelTrackFilter.h"
00019 #include "RecoPixelVertexing/PixelTrackFitting/interface/PixelTrackFilterFactory.h"
00020 
00021 #include "RecoPixelVertexing/PixelTrackFitting/interface/PixelTrackCleaner.h"
00022 #include "RecoPixelVertexing/PixelTrackFitting/interface/PixelTrackCleanerFactory.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 
00029 #include <vector>
00030 
00031 using namespace pixeltrackfitting;
00032 using namespace ctfseeding;
00033 using edm::ParameterSet;
00034 
00035 PixelTrackProducer::PixelTrackProducer(const ParameterSet& cfg)
00036   : theConfig(cfg), theFitter(0), theFilter(0), theCleaner(0), theGenerator(0), theRegionProducer(0)
00037 {
00038   edm::LogInfo("PixelTrackProducer")<<" construction...";
00039   produces<reco::TrackCollection>();
00040   produces<TrackingRecHitCollection>();
00041   produces<reco::TrackExtraCollection>();
00042 }
00043 
00044 
00045 PixelTrackProducer::~PixelTrackProducer()
00046 {
00047 }
00048 
00049 void PixelTrackProducer::endRun(edm::Run &run, const edm::EventSetup& es)
00050 { 
00051   delete theFilter;
00052   delete theFitter;
00053   delete theCleaner;
00054   delete theGenerator;
00055   delete theRegionProducer;
00056 }
00057 
00058 void PixelTrackProducer::beginRun(edm::Run &run, const edm::EventSetup& es)
00059 {
00060   ParameterSet regfactoryPSet = theConfig.getParameter<ParameterSet>("RegionFactoryPSet");
00061   std::string regfactoryName = regfactoryPSet.getParameter<std::string>("ComponentName");
00062   theRegionProducer = TrackingRegionProducerFactory::get()->create(regfactoryName,regfactoryPSet);
00063 
00064   ParameterSet orderedPSet =
00065       theConfig.getParameter<ParameterSet>("OrderedHitsFactoryPSet");
00066   std::string orderedName = orderedPSet.getParameter<std::string>("ComponentName");
00067   theGenerator = OrderedHitsGeneratorFactory::get()->create( orderedName, orderedPSet);
00068 
00069   ParameterSet fitterPSet = theConfig.getParameter<ParameterSet>("FitterPSet");
00070   std::string fitterName = fitterPSet.getParameter<std::string>("ComponentName");
00071   theFitter = PixelFitterFactory::get()->create( fitterName, fitterPSet);
00072 
00073   ParameterSet filterPSet = theConfig.getParameter<ParameterSet>("FilterPSet");
00074   std::string  filterName = filterPSet.getParameter<std::string>("ComponentName");
00075   theFilter = PixelTrackFilterFactory::get()->create( filterName, filterPSet);
00076 
00077   ParameterSet cleanerPSet = theConfig.getParameter<ParameterSet>("CleanerPSet");
00078   std::string  cleanerName = cleanerPSet.getParameter<std::string>("ComponentName");
00079   theCleaner = PixelTrackCleanerFactory::get()->create( cleanerName, cleanerPSet);
00080 
00081 }
00082 
00083 void PixelTrackProducer::produce(edm::Event& ev, const edm::EventSetup& es)
00084 {
00085   LogDebug("PixelTrackProducer, produce")<<"event# :"<<ev.id();
00086 
00087   TracksWithRecHits tracks;
00088 
00089   typedef std::vector<TrackingRegion* > Regions;
00090   typedef Regions::const_iterator IR;
00091   Regions regions = theRegionProducer->regions(ev,es);
00092 
00093   for (IR ir=regions.begin(), irEnd=regions.end(); ir < irEnd; ++ir) {
00094     const TrackingRegion & region = **ir;
00095 
00096     const OrderedSeedingHits & triplets =  theGenerator->run(region,ev,es); 
00097     unsigned int nTriplets = triplets.size();
00098     edm::LogInfo("PixelTrackProducer") << "number of proto tracks: " << triplets.size();
00099 
00100     // producing tracks
00101     for (unsigned int iTriplet = 0; iTriplet < nTriplets; ++iTriplet) { 
00102       const SeedingHitSet & triplet = triplets[iTriplet]; 
00103 
00104       std::vector<const TrackingRecHit *> hits;
00105       for (unsigned int iHit = 0, nHits = triplet.size(); iHit < nHits; ++iHit) {
00106         hits.push_back( triplet[iHit] );
00107       }
00108 
00109       // fitting
00110       reco::Track* track = theFitter->run(es, hits, region);
00111 
00112       // decide if track should be skipped according to filter 
00113       if ( ! (*theFilter)(track) ) { 
00114         delete track; 
00115         continue; 
00116       }
00117 
00118       // add tracks 
00119       tracks.push_back(TrackWithRecHits(track, hits));
00120     }
00121   }
00122 
00123   // skip ovelrapped tracks
00124   if(theCleaner) tracks = theCleaner->cleanTracks(tracks);
00125 
00126   // store tracks
00127   store(ev, tracks);
00128 
00129   // clean memory
00130   for (IR ir=regions.begin(), irEnd=regions.end(); ir < irEnd; ++ir) delete (*ir); 
00131 }
00132 
00133 void PixelTrackProducer::store(edm::Event& ev, const TracksWithRecHits & cleanedTracks)
00134 {
00135   std::auto_ptr<reco::TrackCollection> tracks(new reco::TrackCollection);
00136   std::auto_ptr<TrackingRecHitCollection> recHits(new TrackingRecHitCollection);
00137   std::auto_ptr<reco::TrackExtraCollection> trackExtras(new reco::TrackExtraCollection);
00138   typedef std::vector<const TrackingRecHit *> RecHits;
00139 
00140 
00141   int cc = 0, nTracks = cleanedTracks.size();
00142 
00143   for (int i = 0; i < nTracks; i++)
00144   {
00145     reco::Track* track =  cleanedTracks.at(i).first;
00146     const RecHits & hits = cleanedTracks.at(i).second;
00147 
00148     for (unsigned int k = 0; k < hits.size(); k++)
00149     {
00150       TrackingRecHit *hit = (hits.at(k))->clone();
00151       track->setHitPattern(*hit, k);
00152       recHits->push_back(hit);
00153     }
00154     tracks->push_back(*track);
00155     delete track;
00156 
00157   }
00158 
00159   LogDebug("TrackProducer") << "put the collection of TrackingRecHit in the event" << "\n";
00160   edm::OrphanHandle <TrackingRecHitCollection> ohRH = ev.put( recHits );
00161 
00162 
00163   for (int k = 0; k < nTracks; k++)
00164   {
00165     reco::TrackExtra* theTrackExtra = new reco::TrackExtra();
00166 
00167     //fill the TrackExtra with TrackingRecHitRef
00168     unsigned int nHits = tracks->at(k).numberOfValidHits();
00169     for(unsigned int i = 0; i < nHits; ++i) {
00170       theTrackExtra->add(TrackingRecHitRef(ohRH,cc));
00171       cc++;
00172     }
00173 
00174     trackExtras->push_back(*theTrackExtra);
00175     delete theTrackExtra;
00176   }
00177 
00178   LogDebug("TrackProducer") << "put the collection of TrackExtra in the event" << "\n";
00179   edm::OrphanHandle<reco::TrackExtraCollection> ohTE = ev.put(trackExtras);
00180 
00181   for (int k = 0; k < nTracks; k++)
00182   {
00183     const reco::TrackExtraRef theTrackExtraRef(ohTE,k);
00184     (tracks->at(k)).setExtra(theTrackExtraRef);
00185   }
00186 
00187   ev.put(tracks);
00188 }

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