Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:11:23

0001 #include <memory>
0002 #include "FastSimulation/Tracking/plugins/PixelTracksProducer.h"
0003 
0004 #include "FWCore/Framework/interface/Event.h"
0005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0006 #include "FWCore/Framework/interface/ConsumesCollector.h"
0007 
0008 #include "DataFormats/Common/interface/Handle.h"
0009 #include "DataFormats/Common/interface/OwnVector.h"
0010 
0011 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
0012 
0013 //Pixel Specific stuff
0014 #include "RecoTracker/TkTrackingRegions/interface/TrackingRegionProducer.h"
0015 #include "RecoTracker/TkTrackingRegions/interface/TrackingRegionProducerFactory.h"
0016 
0017 #include "RecoTracker/PixelTrackFitting/interface/PixelFitter.h"
0018 
0019 #include "RecoTracker/PixelTrackFitting/interface/PixelTrackFilter.h"
0020 #include "RecoTracker/PixelTrackFitting/interface/TracksWithHits.h"
0021 #include "RecoTracker/TkTrackingRegions/interface/TrackingRegion.h"
0022 
0023 #include "DataFormats/TrackReco/interface/Track.h"
0024 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0025 #include "DataFormats/TrackReco/interface/TrackExtra.h"
0026 #include "DataFormats/Common/interface/OrphanHandle.h"
0027 #include <vector>
0028 
0029 using namespace pixeltrackfitting;
0030 
0031 PixelTracksProducer::PixelTracksProducer(const edm::ParameterSet& conf)
0032     : theRegionProducer(nullptr), ttopoToken(esConsumes()) {
0033   produces<reco::TrackCollection>();
0034   produces<TrackingRecHitCollection>();
0035   produces<reco::TrackExtraCollection>();
0036 
0037   const edm::ParameterSet& regfactoryPSet = conf.getParameter<edm::ParameterSet>("RegionFactoryPSet");
0038   std::string regfactoryName = regfactoryPSet.getParameter<std::string>("ComponentName");
0039   theRegionProducer = TrackingRegionProducerFactory::get()->create(regfactoryName, regfactoryPSet, consumesCollector());
0040 
0041   fitterToken = consumes<PixelFitter>(conf.getParameter<edm::InputTag>("Fitter"));
0042   filterToken = consumes<PixelTrackFilter>(conf.getParameter<edm::InputTag>("Filter"));
0043 
0044   // The name of the seed producer
0045   auto seedProducer = conf.getParameter<edm::InputTag>("SeedProducer");
0046   seedProducerToken = consumes<TrajectorySeedCollection>(seedProducer);
0047 }
0048 
0049 // Virtual destructor needed.
0050 PixelTracksProducer::~PixelTracksProducer() = default;
0051 
0052 // Functions that gets called by framework every event
0053 void PixelTracksProducer::produce(edm::Event& e, const edm::EventSetup& es) {
0054   std::unique_ptr<reco::TrackCollection> tracks(new reco::TrackCollection);
0055   std::unique_ptr<TrackingRecHitCollection> recHits(new TrackingRecHitCollection);
0056   std::unique_ptr<reco::TrackExtraCollection> trackExtras(new reco::TrackExtraCollection);
0057   typedef std::vector<const TrackingRecHit*> RecHits;
0058 
0059   TracksWithRecHits pixeltracks;
0060   TracksWithRecHits cleanedTracks;
0061 
0062   const TrackerTopology& ttopo = es.getData(ttopoToken);
0063 
0064   edm::Handle<PixelFitter> hfitter;
0065   e.getByToken(fitterToken, hfitter);
0066   const PixelFitter& fitter = *hfitter;
0067 
0068   edm::Handle<PixelTrackFilter> hfilter;
0069   e.getByToken(filterToken, hfilter);
0070   const PixelTrackFilter& theFilter = *hfilter;
0071 
0072   edm::Handle<TrajectorySeedCollection> theSeeds;
0073   e.getByToken(seedProducerToken, theSeeds);
0074 
0075   // No seed -> output an empty track collection
0076   if (theSeeds->empty()) {
0077     e.put(std::move(tracks));
0078     e.put(std::move(recHits));
0079     e.put(std::move(trackExtras));
0080     return;
0081   }
0082 
0083   //only one region Global, but it is called at every event...
0084   //maybe there is a smarter way to set it only once
0085   //NEED TO FIX
0086   typedef std::vector<std::unique_ptr<TrackingRegion> > Regions;
0087   typedef Regions::const_iterator IR;
0088   Regions regions = theRegionProducer->regions(e, es);
0089   for (IR ir = regions.begin(), irEnd = regions.end(); ir < irEnd; ++ir) {
0090     const TrackingRegion& region = **ir;
0091 
0092     // Loop over the seeds
0093     TrajectorySeedCollection::const_iterator aSeed = theSeeds->begin();
0094     TrajectorySeedCollection::const_iterator lastSeed = theSeeds->end();
0095     for (; aSeed != lastSeed; ++aSeed) {
0096       // Loop over the rechits
0097       std::vector<const TrackingRecHit*> TripletHits(3, static_cast<const TrackingRecHit*>(nullptr));
0098       unsigned int iRecHit = 0;
0099       for (auto const& recHit : aSeed->recHits()) {
0100         TripletHits[iRecHit] = &recHit;
0101         ++iRecHit;
0102       }
0103 
0104       // fitting the triplet
0105       std::unique_ptr<reco::Track> track = fitter.run(TripletHits, region);
0106 
0107       // decide if track should be skipped according to filter
0108       if (!theFilter(track.get(), TripletHits)) {
0109         continue;
0110       }
0111 
0112       // add tracks
0113       pixeltracks.push_back(TrackWithRecHits(track.release(), TripletHits));
0114     }
0115   }
0116 
0117   int cc = 0;
0118   int nTracks = pixeltracks.size();
0119   for (int i = 0; i < nTracks; ++i) {
0120     reco::Track* track = pixeltracks.at(i).first;
0121     const RecHits& hits = pixeltracks.at(i).second;
0122 
0123     for (unsigned int k = 0; k < hits.size(); k++) {
0124       TrackingRecHit* hit = (hits.at(k))->clone();
0125       track->appendHitPattern(*hit, ttopo);
0126       recHits->push_back(hit);
0127     }
0128 
0129     tracks->push_back(*track);
0130     delete track;
0131   }
0132 
0133   edm::OrphanHandle<TrackingRecHitCollection> ohRH = e.put(std::move(recHits));
0134   edm::RefProd<TrackingRecHitCollection> ohRHProd(ohRH);
0135 
0136   for (int k = 0; k < nTracks; ++k) {
0137     // reco::TrackExtra* theTrackExtra = new reco::TrackExtra();
0138     reco::TrackExtra theTrackExtra;
0139 
0140     //fill the TrackExtra with TrackingRecHitRef
0141     // unsigned int nHits = tracks->at(k).numberOfValidHits();
0142     const unsigned nHits = 3;  // We are dealing with triplets!
0143     theTrackExtra.setHits(ohRHProd, cc, nHits);
0144     cc += nHits;
0145 
0146     trackExtras->push_back(theTrackExtra);
0147     //trackExtras->push_back(*theTrackExtra);
0148     //delete theTrackExtra;
0149   }
0150 
0151   edm::OrphanHandle<reco::TrackExtraCollection> ohTE = e.put(std::move(trackExtras));
0152 
0153   for (int k = 0; k < nTracks; k++) {
0154     const reco::TrackExtraRef theTrackExtraRef(ohTE, k);
0155     (tracks->at(k)).setExtra(theTrackExtraRef);
0156   }
0157 
0158   e.put(std::move(tracks));
0159 }
0160 
0161 #include "FWCore/Framework/interface/MakerMacros.h"
0162 DEFINE_FWK_MODULE(PixelTracksProducer);