Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:28:39

0001 #include <vector>
0002 
0003 #include "DataFormats/TrackReco/interface/Track.h"
0004 #include "DataFormats/TrackReco/interface/TrackExtra.h"
0005 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0006 #include "FWCore/Framework/interface/Event.h"
0007 #include "FWCore/Framework/interface/EventSetup.h"
0008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0009 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0011 #include "RecoTracker/PixelTrackFitting/interface/PixelFitter.h"
0012 #include "RecoTracker/PixelTrackFitting/interface/PixelTrackCleanerWrapper.h"
0013 #include "RecoTracker/PixelTrackFitting/interface/PixelTrackFilter.h"
0014 #include "RecoTracker/PixelTrackFitting/interface/PixelTrackReconstruction.h"
0015 #include "RecoTracker/TkHitPairs/interface/RegionsSeedingHitSets.h"
0016 #include "RecoTracker/TkTrackingRegions/interface/TrackingRegion.h"
0017 
0018 using namespace pixeltrackfitting;
0019 using edm::ParameterSet;
0020 
0021 PixelTrackReconstruction::PixelTrackReconstruction(const ParameterSet& cfg, edm::ConsumesCollector&& iC)
0022     : theHitSetsToken(iC.consumes<RegionsSeedingHitSets>(cfg.getParameter<edm::InputTag>("SeedingHitSets"))),
0023       theFitterToken(iC.consumes<PixelFitter>(cfg.getParameter<edm::InputTag>("Fitter"))) {
0024   edm::InputTag filterTag = cfg.getParameter<edm::InputTag>("Filter");
0025   if (not filterTag.label().empty()) {
0026     theFilterToken = iC.consumes<PixelTrackFilter>(filterTag);
0027   }
0028   std::string cleanerName = cfg.getParameter<std::string>("Cleaner");
0029   if (not cleanerName.empty()) {
0030     theCleanerToken = iC.esConsumes(edm::ESInputTag("", cleanerName));
0031   }
0032 }
0033 
0034 PixelTrackReconstruction::~PixelTrackReconstruction() {}
0035 
0036 void PixelTrackReconstruction::fillDescriptions(edm::ParameterSetDescription& desc) {
0037   desc.add<edm::InputTag>("SeedingHitSets", edm::InputTag("pixelTracksHitTriplets"));
0038   desc.add<edm::InputTag>("Fitter", edm::InputTag("pixelFitterByHelixProjections"));
0039   desc.add<edm::InputTag>("Filter", edm::InputTag("pixelTrackFilterByKinematics"));
0040   desc.add<std::string>("Cleaner", "pixelTrackCleanerBySharedHits");
0041 }
0042 
0043 void PixelTrackReconstruction::run(TracksWithTTRHs& tracks, edm::Event& ev, const edm::EventSetup& es) {
0044   edm::Handle<RegionsSeedingHitSets> hhitSets;
0045   ev.getByToken(theHitSetsToken, hhitSets);
0046   const auto& hitSets = *hhitSets;
0047 
0048   edm::Handle<PixelFitter> hfitter;
0049   ev.getByToken(theFitterToken, hfitter);
0050   const auto& fitter = *hfitter;
0051 
0052   const PixelTrackFilter* filter = nullptr;
0053   if (!theFilterToken.isUninitialized()) {
0054     edm::Handle<PixelTrackFilter> hfilter;
0055     ev.getByToken(theFilterToken, hfilter);
0056     filter = hfilter.product();
0057   }
0058 
0059   std::vector<const TrackingRecHit*> hits;
0060   hits.reserve(4);
0061   for (const auto& regionHitSets : hitSets) {
0062     const TrackingRegion& region = regionHitSets.region();
0063 
0064     for (const SeedingHitSet& tuplet : regionHitSets) {
0065       /// FIXME at some point we need to migrate the fitter...
0066       auto nHits = tuplet.size();
0067       hits.resize(nHits);
0068       for (unsigned int iHit = 0; iHit < nHits; ++iHit)
0069         hits[iHit] = tuplet[iHit];
0070 
0071       // fitting
0072       std::unique_ptr<reco::Track> track = fitter.run(hits, region);
0073       if (!track)
0074         continue;
0075 
0076       if (filter) {
0077         if (!(*filter)(track.get(), hits)) {
0078           continue;
0079         }
0080       }
0081       /* roll back to verify if HLT tau is affected  
0082       // all legacy tracks are "highPurity"
0083       // setting all others bits as well (as in ckf)
0084       track->setQuality(reco::TrackBase::loose);
0085       track->setQuality(reco::TrackBase::tight);
0086       track->setQuality(reco::TrackBase::highPurity);
0087       */
0088       // add tracks
0089       tracks.emplace_back(track.release(), tuplet);
0090     }
0091   }
0092 
0093   // skip ovelrapped tracks
0094   if (theCleanerToken.isInitialized()) {
0095     const auto& cleaner = es.getData(theCleanerToken);
0096     if (cleaner.fast())
0097       cleaner.cleanTracks(tracks);
0098     else
0099       tracks = PixelTrackCleanerWrapper(&cleaner).clean(tracks);
0100   }
0101 }